setTurbulenceFields.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) 2022 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26Application
27 setTurbulenceFields
28
29Group
30 grpPreProcessingUtilities
31
32Description
33 Initialises turbulence fields according to
34 various empirical governing equations.
35
36 References:
37 \verbatim
38 Initialisation method (tag:M):
39 Manceau, R. (n.d.).
40 A two-step automatic initialization
41 procedure for RANS computations.
42 (Unpublished).
43
44 Previous-version of the initialisation model (tag:LM):
45 Lardeau, S., & Manceau, R. (2014).
46 Computations of complex flow configurations using
47 a modified elliptic-blending Reynolds-stress model.
48 10th International ERCOFTAC Symposium on Engineering
49 Turbulence Modelling and Measurements. Marbella, Spain.
50 https://hal.archives-ouvertes.fr/hal-01051799
51 \endverbatim
52
53Usage
54 Minimal example by using \c system/setTurbulenceFieldsDict:
55 \verbatim
56 // Mandatory entries
57 uRef <scalar>;
58
59 // Optional entries
60 initialiseU <bool>;
61 initialiseEpsilon <bool>;
62 initialiseK <bool>;
63 initialiseOmega <bool>;
64 initialiseR <bool>;
65 writeF <bool>;
66
67 kappa <scalar>;
68 Cmu <scalar>;
69 dPlusRef <scalar>;
70
71 f <word>;
72 U <word>;
73 epsilon <word>;
74 k <word>;
75 omega <word>;
76 R <word>;
77 \endverbatim
78
79 where the entries mean:
80 \table
81 Property | Description | Type | Reqd | Deflt
82 uRef | Reference speed | scalar | yes | -
83 initialiseU | Flag to initialise U | bool | no | false
84 initialiseEpsilon | Flag to initialise epsilon | bool | no | false
85 initialiseK | Flag to initialise k | bool | no | false
86 initialiseOmega | Flag to initialise omega | bool | no | false
87 initialiseR | Flag to initialise R | bool | no | false
88 writeF | Flag to write elliptic-blending field, f | bool | no | false
89 kappa | von Karman constant | scalar | no | 0.41
90 Cmu | Empirical constant | scalar | no | 0.09
91 dPlusRef | Reference dPlus | scalar | no | 15
92 f | Name of operand f field | word | no | f
93 U | Name of operand U field | word | no | U
94 epsilon | Name of operand epsilon field | word | no | epsilon
95 k | Name of operand k field | word | no | k
96 omega | Name of operand omega field | word | no | omega
97 R | Name of operand R field | word | no | R
98 \endtable
99
100Note
101 - Time that the utility applies to is determined by the
102 \c controlDict.startFrom and \c controlDict.startTime entries.
103 - The utility modifies near-wall fields, hence
104 can be more effective for low-Re mesh cases.
105
106\*---------------------------------------------------------------------------*/
107
108#include "fvCFD.H"
112#include "wallFvPatch.H"
114
115using namespace Foam;
116
117// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
118
119void InfoField(const word& fldName)
120{
121 Info<< "Writing field: " << fldName << nl << endl;
122}
123
124
125IOobject createIOobject
126(
127 const fvMesh& mesh,
128 const word& name,
129 IOobject::readOption rOpt = IOobject::READ_IF_PRESENT
130)
131{
132 return
134 (
135 name,
136 mesh.time().timeName(),
137 mesh,
138 rOpt,
139 IOobject::NO_WRITE,
140 false // do not register
141 );
142}
143
144
146(
147 const fvMesh& mesh,
148 const volVectorField& U
149)
150{
151 const Time& runTime = mesh.time();
152
153 if
154 (
156 (
157 basicThermo::dictName,
158 runTime.constant(),
159 mesh
161 )
162 {
163 // Compressible
164 autoPtr<fluidThermo> pThermo(fluidThermo::New(mesh));
167
168 // Create phi
169 #include "compressibleCreatePhi.H"
170
172 (
173 compressible::turbulenceModel::New
174 (
175 rho,
176 U,
177 phi,
178 thermo
179 )
180 );
181
182 turbulence->validate();
183
185 }
186 else
187 {
188 // Incompressible
189
190 // Create phi
191 #include "createPhi.H"
192
194
196 (
197 incompressible::turbulenceModel::New(U, phi, laminarTransport)
198 );
199
200 turbulence->validate();
201
203 }
204}
205
206
207// Calculate elliptic blending function
208// between near-wall and weakly-inhomogeneous regions
209void calcF
210(
211 const volScalarField& L,
213)
214{
215 tmp<volScalarField> tinvLsqr = scalar(1)/sqr(L);
216 const volScalarField& invLsqr = tinvLsqr.cref();
217
218 // (M:Eq. 6)
220 (
221 fvm::Sp(invLsqr, f)
222 - fvm::laplacian(f)
223 ==
224 invLsqr
225 );
226
227 tinvLsqr.clear();
228
229 fEqn.ref().relax();
230 solve(fEqn);
231
232 // (M:p. 2)
233 const dimensioned<scalarMinMax> fMinMax
234 (
235 dimless,
236 scalarMinMax(Zero, scalar(1) - Foam::exp(-scalar(400)/scalar(50)))
237 );
238
239 f.clip(fMinMax);
240}
241
242
243// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
244// Main program:
245int main(int argc, char *argv[])
246{
247 argList::addNote
248 (
249 "Sets initial turbulence fields based on"
250 " various empirical equations"
251 );
252 argList::noFunctionObjects();
253 argList::addOption("dict", "file", "Alternative setTurbulenceFieldsDict");
254
255 #include "addRegionOption.H"
256
257 #include "setRootCase.H"
258 #include "createTime.H"
259
260 const word dictName("setTurbulenceFieldsDict");
262 Info<< "Reading " << dictIO.name() << nl << endl;
264
265 #include "createNamedMesh.H"
266
267 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
268
269 IOstream::defaultPrecision(15);
270
271
272 // Dictionary input (M:p. 2)
273
274 const scalar uRef = dict.getCheck<scalar>("uRef", scalarMinMax::ge(SMALL));
275 const dimensionedScalar uTau(dimVelocity, 0.05*uRef);
276 const scalar kappa =
277 dict.getCheckOrDefault<scalar>("kappa", 0.41, scalarMinMax::ge(SMALL));
278 const scalar Cmu =
279 dict.getCheckOrDefault<scalar>("Cmu", 0.09, scalarMinMax::ge(SMALL));
280 const scalar dPlusRef =
281 dict.getCheckOrDefault<scalar>("dPlusRef", 15, scalarMinMax::ge(SMALL));
282 const word fName = dict.getOrDefault<word>("f", "f");
283 const word UName = dict.getOrDefault<word>("U", "U");
284 const word epsilonName = dict.getOrDefault<word>("epsilon", "epsilon");
285 const word kName = dict.getOrDefault<word>("k", "k");
286 const word omegaName = dict.getOrDefault<word>("omega", "omega");
287 const word RName = dict.getOrDefault<word>("R", "R");
288 const bool initU = dict.getOrDefault<bool>("initialiseU", false);
289 const bool initEpsilon =
290 dict.getOrDefault<bool>("initialiseEpsilon", false);
291 const bool initK = dict.getOrDefault<bool>("initialiseK", false);
292 const bool initOmega = dict.getOrDefault<bool>("initialiseOmega", false);
293 const bool initR = dict.getOrDefault<bool>("initialiseR", false);
294 const bool writeF = dict.getOrDefault<bool>("writeF", false);
295
296
297 // Start initialising the operand fields
298
299 // Read operand fields
300
302 (
303 createIOobject(mesh, UName, IOobject::MUST_READ),
304 mesh
305 );
306
307 // Infer the initial BCs from the velocity
308 const wordList bcTypes
309 (
310 tU.cref().boundaryField().size(),
312 );
313
314 tmp<volScalarField> tepsilon;
316 tmp<volScalarField> tomega;
318
319 if (initEpsilon)
320 {
321 tepsilon = tmp<volScalarField>::New
322 (
323 createIOobject(mesh, epsilonName),
324 mesh,
325 dimensionedScalar(sqr(dimLength)/pow3(dimTime), SMALL),
326 bcTypes
327 );
328 }
329
330 if (initK)
331 {
333 (
334 createIOobject(mesh, kName),
335 mesh,
336 dimensionedScalar(sqr(dimLength/dimTime), SMALL),
337 bcTypes
338 );
339 }
340
341 if (initOmega)
342 {
344 (
345 createIOobject(mesh, omegaName),
346 mesh,
347 dimensionedScalar(dimless/dimTime, SMALL),
348 bcTypes
349 );
350 }
351
352 if (initR)
353 {
355 (
356 createIOobject(mesh, RName),
357 mesh,
358 dimensionedSymmTensor(sqr(dimLength/dimTime), Zero),
359 bcTypes
360 );
361 }
362
363
364 // Create elliptic blending factor field
365
367 (
369 (
370 fName,
371 runTime.timeName(),
372 mesh,
373 IOobject::NO_READ,
374 IOobject::NO_WRITE,
375 false // do not register
376 ),
377 mesh,
378 dimensionedScalar(dimless, scalar(1)),
380 );
381
382 for (fvPatchScalarField& pfld : f.boundaryFieldRef())
383 {
384 if (isA<wallFvPatch>(pfld.patch()))
385 {
386 pfld == Zero;
387 }
388 }
389
390
391 // Create auxillary fields for the initialisation
392
393 tmp<volScalarField> tnu = nu(mesh, tU.cref());
394 const volScalarField& nu = tnu.cref();
395
396 // (M:p. 2)
397 tmp<volScalarField> tL = scalar(50)*nu/uTau;
398 const volScalarField& L = tL.cref();
399
400 calcF(L, f);
401
402 // (M:Eq. 8)
403 tmp<volScalarField> td = -tL*Foam::log(scalar(1) - f);
404 const volScalarField& d = td.cref();
405
406 // (M:p. 2)
407 const volScalarField dPlus(d*uTau/nu);
408
409 // (M:Eq. 11)
411 (
412 pow4(uTau)/(kappa*nu*max(dPlus, dPlusRef))
413 );
414
415 // (M:Eq. 13)
416 const volScalarField fk(Foam::exp(-dPlus/scalar(25)));
417
418 // (M:Eq. 12)
419 const volScalarField k
420 (
421 (epsilon*sqr(td)*sqr(fk))/(2*nu)
422 + sqr(uTau)*sqr(scalar(1) - fk)/Foam::sqrt(Cmu)
423 );
424
425
426 // Initialise operand fields
427
428 if (initU)
429 {
430 volVectorField& U = tU.ref();
431
432 // Reichardt’s law (M:Eq. 10)
433 const scalar C = 7.8;
434 const scalar B1 = 11;
435 const scalar B2 = 3;
436 const volScalarField fRei
437 (
438 Foam::log(scalar(1) + kappa*dPlus)/kappa
439 + C*
440 (
441 scalar(1)
442 - Foam::exp(-dPlus/B1)
443 - dPlus/B1*Foam::exp(-dPlus/B2)
444 )
445 );
446
447 // (M:Eq. 9)
448 const dimensionedScalar maxU(dimVelocity, SMALL);
449 U *= min(scalar(1), fRei*uTau/max(mag(U), maxU));
450 }
451
452 if (tepsilon.valid())
453 {
454 tepsilon.ref() = epsilon;
455 }
456
457 if (tk.valid())
458 {
459 tk.ref() = k;
460 }
461
462 if (tomega.valid())
463 {
464 const dimensionedScalar k0(sqr(dimLength/dimTime), SMALL);
465 tomega.ref() = Cmu*epsilon/(k + k0);
466 }
467
468 if (tR.valid())
469 {
470 volSymmTensorField& R = tR.ref();
471
472 // (M:Eq. 3)
473 const volSphericalTensorField Rdiag(k*twoThirdsI);
474 forAll(R, celli)
475 {
476 R[celli] = Rdiag[celli];
477 }
478 }
479
480
481 // Write operand fields
482
483 Info<< endl;
484
485 if (initU)
486 {
487 InfoField(tU->name());
488 tU->write();
489 }
490
491 if (tepsilon.valid())
492 {
493 InfoField(tepsilon->name());
494 tepsilon->write();
495 }
496
497 if (tk.valid())
498 {
499 InfoField(tk->name());
500 tk->write();
501 }
502
503 if (tomega.valid())
504 {
505 InfoField(tomega->name());
506 tomega->write();
507 }
508
509 if (tR.valid())
510 {
511 InfoField(tR->name());
512 tR->write();
513 }
514
515 if (writeF)
516 {
517 InfoField(f.name());
518 f.write();
519 }
520
521
522 Info<< nl;
523 runTime.printExecutionTime(Info);
524
525 Info<< "End\n" << endl;
526
527 return 0;
528}
529
530
531// ************************************************************************* //
label k
#define R(A, B, C, D, E, F, K, M)
Y[inertIndex] max(0.0)
surfaceScalarField & phi
Graphite solid properties.
Definition: C.H:53
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
bool typeHeaderOk(const bool checkType=true, const bool search=true, const bool verbose=true)
Read header (uses typeFilePath to find file) and check its info.
readOption
Enumeration defining the read options.
Definition: IOobject.H:177
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
Generic dimensioned Type class.
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
Fundamental fluid thermodynamic properties.
Definition: fluidThermo.H:56
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
A simple single-phase transport model based on viscosityModel.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
A class for managing temporary objects.
Definition: tmp.H:65
const T & cref() const
Definition: tmpI.H:213
void clear() const noexcept
Definition: tmpI.H:287
bool valid() const noexcept
Identical to good(), or bool operator.
Definition: tmp.H:292
T & ref() const
Definition: tmpI.H:227
A class for handling words, derived from Foam::string.
Definition: word.H:68
U
Definition: pEqn.H:72
Creates and initialises the face-flux field phi.
dynamicFvMesh & mesh
engineTime & runTime
Required Variables.
scalar uTau
scalar epsilon
const word dictName("faMeshDefinition")
compressible::turbulenceModel & turbulence
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
Namespace for OpenFOAM.
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
messageStream Info
Information stream (stdout output on master, null elsewhere)
dimensionedScalar log(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
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
dimensionedScalar pow4(const dimensionedScalar &ds)
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
labelList f(nPoints)
volScalarField & nu
dictionary dict
IOobject dictIO
CEqn solve()
Info<< "Reading thermophysical properties\n"<< endl;autoPtr< psiReactionThermo > pThermo(psiReactionThermo::New(mesh))
singlePhaseTransportModel laminarTransport(U, phi)
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
static const char *const typeName
The type name used in ensight case files.
const word UName(propsDict.getOrDefault< word >("U", "U"))
const vector L(dict.get< vector >("L"))