Loading...
Searching...
No Matches
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-2023 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"
113#include "processorFvPatch.H"
115
116using namespace Foam;
117
118// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
119
120void InfoField(const word& fldName)
121{
122 Info<< "Writing field: " << fldName << nl << endl;
123}
124
125
126template<class Type>
127void correctProcessorPatches(GeometricField<Type, fvPatchField, volMesh>& fld)
128{
129 if (UPstream::parRun())
130 {
131 fld.boundaryFieldRef().template evaluateCoupled<processorFvPatch>();
132 }
133}
134
135
136IOobject createIOobject
137(
138 const fvMesh& mesh,
139 const word& name,
141)
142{
143 return
145 (
146 name,
147 mesh.time().timeName(),
148 mesh,
149 rOpt,
152 );
153}
154
155
157(
158 const fvMesh& mesh,
159 const volVectorField& U
160)
161{
162 const Time& runTime = mesh.time();
163
164 if
165 (
167 (
169 runTime.constant(),
170 mesh
172 )
173 {
174 // Compressible
178
179 // Create phi
180 #include "compressibleCreatePhi.H"
181
183 (
185 (
186 rho,
187 U,
188 phi,
189 thermo
190 )
191 );
192
193 turbulence->validate();
194
196 }
197 else
198 {
199 // Incompressible
200
201 // Create phi
202 #include "createPhi.H"
203
205
207 (
209 );
210
211 turbulence->validate();
212
214 }
215}
216
217
218// Calculate elliptic blending function
219// between near-wall and weakly-inhomogeneous regions
220void calcF
221(
222 const volScalarField& L,
224)
225{
226 tmp<volScalarField> tinvLsqr = scalar(1)/sqr(L);
227 const volScalarField& invLsqr = tinvLsqr.cref();
228
229 // (M:Eq. 6)
231 (
232 fvm::Sp(invLsqr, f)
234 ==
235 invLsqr
236 );
237
238 tinvLsqr.clear();
239
240 fEqn.ref().relax();
241 solve(fEqn);
242
243 // (M:p. 2)
244 f.clamp_range(0, scalar(1) - Foam::exp(-scalar(400)/scalar(50)));
245}
246
247
248// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
249// Main program:
250int main(int argc, char *argv[])
251{
253 (
254 "Sets initial turbulence fields based on"
255 " various empirical equations"
256 );
258 argList::addOption("dict", "file", "Alternative setTurbulenceFieldsDict");
259
260 #include "addRegionOption.H"
261
262 #include "setRootCase.H"
263 #include "createTime.H"
264
265 const word dictName("setTurbulenceFieldsDict");
267 Info<< "Reading " << dictIO.name() << nl << endl;
269
270 #include "createNamedMesh.H"
271
272 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
273
275
276
277 // Dictionary input (M:p. 2)
278
279 const scalar uRef = dict.getCheck<scalar>("uRef", scalarMinMax::ge(SMALL));
280 const dimensionedScalar uTau(dimVelocity, 0.05*uRef);
281 const scalar kappa =
282 dict.getCheckOrDefault<scalar>("kappa", 0.41, scalarMinMax::ge(SMALL));
283 const scalar Cmu =
284 dict.getCheckOrDefault<scalar>("Cmu", 0.09, scalarMinMax::ge(SMALL));
285 const scalar dPlusRef =
286 dict.getCheckOrDefault<scalar>("dPlusRef", 15, scalarMinMax::ge(SMALL));
287 const word fName = dict.getOrDefault<word>("f", "f");
288 const word UName = dict.getOrDefault<word>("U", "U");
289 const word epsilonName = dict.getOrDefault<word>("epsilon", "epsilon");
290 const word kName = dict.getOrDefault<word>("k", "k");
291 const word omegaName = dict.getOrDefault<word>("omega", "omega");
292 const word RName = dict.getOrDefault<word>("R", "R");
293 const bool initU = dict.getOrDefault<bool>("initialiseU", false);
294 const bool initEpsilon =
295 dict.getOrDefault<bool>("initialiseEpsilon", false);
296 const bool initK = dict.getOrDefault<bool>("initialiseK", false);
297 const bool initOmega = dict.getOrDefault<bool>("initialiseOmega", false);
298 const bool initR = dict.getOrDefault<bool>("initialiseR", false);
299 const bool writeF = dict.getOrDefault<bool>("writeF", false);
300
301
302 // Start initialising the operand fields
303
304 // Read operand fields
305
307 (
308 createIOobject(mesh, UName, IOobject::MUST_READ),
309 mesh
310 );
311
312 // Infer the initial BCs from the velocity
313 const wordList bcTypes
314 (
315 tU.cref().boundaryField().size(),
316 fixedValueFvPatchScalarField::typeName
317 );
318
319 tmp<volScalarField> tepsilon;
321 tmp<volScalarField> tomega;
323
324 if (initEpsilon)
325 {
326 tepsilon = tmp<volScalarField>::New
327 (
328 createIOobject(mesh, epsilonName),
329 mesh,
331 bcTypes
332 );
333 }
334
335 if (initK)
336 {
338 (
339 createIOobject(mesh, kName),
340 mesh,
342 bcTypes
343 );
344 }
345
346 if (initOmega)
347 {
349 (
350 createIOobject(mesh, omegaName),
351 mesh,
353 bcTypes
354 );
355 }
356
357 if (initR)
358 {
360 (
361 createIOobject(mesh, RName),
362 mesh,
364 bcTypes
365 );
366 }
367
368
369 // Create elliptic blending factor field
370
372 (
374 (
375 fName,
376 runTime.timeName(),
377 mesh,
381 ),
382 mesh,
383 dimensionedScalar(dimless, scalar(1)),
384 fixedValueFvPatchScalarField::typeName
385 );
386
387 for (fvPatchScalarField& pfld : f.boundaryFieldRef())
388 {
389 if (isA<wallFvPatch>(pfld.patch()))
390 {
391 pfld == Zero;
392 }
393 }
394
395
396 // Create auxillary fields for the initialisation
397
398 tmp<volScalarField> tnu = nu(mesh, tU.cref());
399 const volScalarField& nu = tnu.cref();
400
401 // (M:p. 2)
402 tmp<volScalarField> tL = scalar(50)*nu/uTau;
403 const volScalarField& L = tL.cref();
404
405 calcF(L, f);
406
407 // (M:Eq. 8)
408 tmp<volScalarField> td = -tL*Foam::log(scalar(1) - f);
409 const volScalarField& d = td.cref();
410
411 // (M:p. 2)
412 const volScalarField dPlus(d*uTau/nu);
413
414 // (M:Eq. 11)
416 (
417 pow4(uTau)/(kappa*nu*max(dPlus, dPlusRef))
418 );
419
420 // (M:Eq. 13)
421 const volScalarField fk(Foam::exp(-dPlus/scalar(25)));
422
423 // (M:Eq. 12)
424 const volScalarField k
425 (
426 (epsilon*sqr(td)*sqr(fk))/(2*nu)
427 + sqr(uTau)*sqr(scalar(1) - fk)/Foam::sqrt(Cmu)
428 );
429
430
431 // Initialise operand fields
432
433 if (initU)
434 {
435 volVectorField& U = tU.ref();
436
437 // Reichardt’s law (M:Eq. 10)
438 const scalar C = 7.8;
439 const scalar B1 = 11;
440 const scalar B2 = 3;
441 const volScalarField fRei
442 (
443 Foam::log(scalar(1) + kappa*dPlus)/kappa
444 + C*
445 (
446 scalar(1)
447 - Foam::exp(-dPlus/B1)
448 - dPlus/B1*Foam::exp(-dPlus/B2)
449 )
450 );
451
452 // (M:Eq. 9)
453 const dimensionedScalar maxU(dimVelocity, SMALL);
454 U *= min(scalar(1), fRei*uTau/max(mag(U), maxU));
455 correctProcessorPatches(U);
456 }
457
458 if (tepsilon.valid())
459 {
460 tepsilon.ref() = epsilon;
461 correctProcessorPatches(tepsilon.ref());
462 }
463
464 if (tk.valid())
465 {
466 tk.ref() = k;
467 correctProcessorPatches(tk.ref());
468 }
469
470 if (tomega.valid())
471 {
472 const dimensionedScalar k0(sqr(dimLength/dimTime), SMALL);
473 tomega.ref() = Cmu*epsilon/(k + k0);
474 correctProcessorPatches(tomega.ref());
475 }
476
477 if (tR.valid())
478 {
479 auto& R = tR.ref();
480
481 // (M:Eq. 3)
483 forAll(R, celli)
484 {
485 R[celli] = Rdiag[celli];
486 }
487 correctProcessorPatches(R);
488 }
489
490
491 // Write operand fields
492
493 Info<< endl;
494
495 if (initU)
496 {
497 InfoField(tU->name());
498 tU->write();
499 }
500
501 if (tepsilon.valid())
502 {
503 InfoField(tepsilon->name());
504 tepsilon->write();
505 }
506
507 if (tk.valid())
508 {
509 InfoField(tk->name());
510 tk->write();
511 }
512
513 if (tomega.valid())
514 {
515 InfoField(tomega->name());
516 tomega->write();
517 }
518
519 if (tR.valid())
520 {
521 InfoField(tR->name());
522 tR->write();
523 }
524
525 if (writeF)
526 {
527 InfoField(f.name());
528 f.write();
529 }
530
531
532 Info<< nl;
533 runTime.printExecutionTime(Info);
534
535 Info<< "End\n" << endl;
536
537 return 0;
538}
539
540
541// ************************************************************************* //
static const Foam::dimensionedScalar C("", Foam::dimTemperature, 234.5)
label k
#define R(A, B, C, D, E, F, K, M)
Required Classes.
const word UName(propsDict.getOrDefault< word >("U", "U"))
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
Generic GeometricField class.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
@ NO_REGISTER
Do not request registration (bool: false).
readOption
Enumeration defining read preferences.
@ NO_READ
Nothing to be read.
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ MUST_READ
Reading required.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
bool typeHeaderOk(const bool checkType=true, const bool search=true, const bool verbose=true)
Read header (respects is_globalIOobject trait) and check its info. A void type suppresses trait and t...
static unsigned int defaultPrecision() noexcept
Return the default precision.
Definition IOstream.H:437
static autoPtr< IncompressibleTurbulenceModel > New(const volVectorField &U, const surfaceScalarField &phi, const transportModel &transportModel, const word &propertiesName=turbulenceModel::propertiesName)
static MinMax< scalar > ge(const scalar &minVal)
static autoPtr< ThermalDiffusivity > New(const alphaField &alpha, const volScalarField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transportModel, const word &propertiesName=turbulenceModel::propertiesName)
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
static bool parRun(const bool on) noexcept
Set as parallel run on/off.
Definition UPstream.H:1669
static void noFunctionObjects(bool addWithOption=false)
Remove '-noFunctionObjects' option and ignore any occurrences.
Definition argList.C:562
static void addOption(const word &optName, const string &param="", const string &usage="", bool advanced=false)
Add an option to validOptions with usage information.
Definition argList.C:400
static void addNote(const string &note)
Add extra notes for the usage information.
Definition argList.C:477
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
static const word dictName
The dictionary name ("thermophysicalProperties").
Fundamental fluid thermodynamic properties.
Definition fluidThermo.H:52
static autoPtr< fluidThermo > New(const fvMesh &, const word &phaseName=word::null)
Selector.
Definition fluidThermo.C:68
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
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:75
const T & cref() const
Return const reference to the object or to the contents of a (non-null) managed pointer.
Definition tmpI.H:221
void clear() const noexcept
If object pointer points to valid object: delete object and set pointer to nullptr.
Definition tmpI.H:289
bool valid() const noexcept
Identical to good(), or bool operator.
Definition tmp.H:481
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
A class for handling words, derived from Foam::string.
Definition word.H:66
U
Definition pEqn.H:72
Creates and initialises the face-flux field phi.
dynamicFvMesh & mesh
engineTime & runTime
Required Classes.
scalar uTau
scalar epsilon
const word dictName("faMeshDefinition")
compressible::turbulenceModel & turbulence
auto & name
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
Namespace for OpenFOAM.
List< word > wordList
List of word.
Definition fileName.H:60
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< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const dimensionSet dimVelocity
messageStream Info
Information stream (stdout output on master, null elsewhere).
dimensionedScalar log(const dimensionedScalar &ds)
dimensioned< symmTensor > dimensionedSymmTensor
Dimensioned tensor obtained from generic dimensioned type.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensionedScalar sqrt(const dimensionedScalar &ds)
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
GeometricField< sphericalTensor, fvPatchField, volMesh > volSphericalTensorField
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
dimensionedScalar pow4(const dimensionedScalar &ds)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
static const sphericalTensor twoThirdsI(2.0/3.0)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
fvPatchField< scalar > fvPatchScalarField
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
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:299
const vector L(dict.get< vector >("L"))