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 -------------------------------------------------------------------------------
10 License
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 
26 Application
27  setTurbulenceFields
28 
29 Group
30  grpPreProcessingUtilities
31 
32 Description
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 
53 Usage
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 
100 Note
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"
110 #include "turbulentTransportModel.H"
112 #include "processorFvPatchField.H"
113 #include "wallFvPatch.H"
114 #include "fixedValueFvPatchFields.H"
115 
116 using namespace Foam;
117 
118 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
119 
120 void InfoField(const word& fldName)
121 {
122  Info<< "Writing field: " << fldName << nl << endl;
123 }
124 
125 
126 template<class Type>
127 void correctProcessorPatches
128 (
130 )
131 {
132  if (!Pstream::parRun())
133  {
134  return;
135  }
136 
137  // Not possible to use correctBoundaryConditions on fields as they may
138  // use local info as opposed to the constraint values employed here,
139  // but still need to update processor patches
140  auto& bf = vf.boundaryFieldRef();
141 
142  forAll(bf, patchi)
143  {
144  if (isA<processorFvPatchField<Type>>(bf[patchi]))
145  {
146  bf[patchi].initEvaluate();
147  }
148  }
149 
150  forAll(bf, patchi)
151  {
152  if (isA<processorFvPatchField<Type>>(bf[patchi]))
153  {
154  bf[patchi].evaluate();
155  }
156  }
157 }
158 
159 
160 IOobject createIOobject
161 (
162  const fvMesh& mesh,
163  const word& name,
165 )
166 {
167  return
168  IOobject
169  (
170  name,
171  mesh.time().timeName(),
172  mesh,
173  rOpt,
175  false // do not register
176  );
177 }
178 
179 
181 (
182  const fvMesh& mesh,
183  const volVectorField& U
184 )
185 {
186  const Time& runTime = mesh.time();
187 
188  if
189  (
190  IOobject
191  (
193  runTime.constant(),
194  mesh
195  ).typeHeaderOk<IOdictionary>(true)
196  )
197  {
198  // Compressible
201  volScalarField rho(thermo.rho());
202 
203  // Create phi
204  #include "compressibleCreatePhi.H"
205 
207  (
209  (
210  rho,
211  U,
212  phi,
213  thermo
214  )
215  );
216 
217  turbulence->validate();
218 
219  return tmp<volScalarField>::New(turbulence->nu());
220  }
221  else
222  {
223  // Incompressible
224 
225  // Create phi
226  #include "createPhi.H"
227 
229 
231  (
233  );
234 
235  turbulence->validate();
236 
237  return tmp<volScalarField>::New(turbulence->nu());
238  }
239 }
240 
241 
242 // Calculate elliptic blending function
243 // between near-wall and weakly-inhomogeneous regions
244 void calcF
245 (
246  const volScalarField& L,
248 )
249 {
250  tmp<volScalarField> tinvLsqr = scalar(1)/sqr(L);
251  const volScalarField& invLsqr = tinvLsqr.cref();
252 
253  // (M:Eq. 6)
255  (
256  fvm::Sp(invLsqr, f)
257  - fvm::laplacian(f)
258  ==
259  invLsqr
260  );
261 
262  tinvLsqr.clear();
263 
264  fEqn.ref().relax();
265  solve(fEqn);
266 
267  // (M:p. 2)
268  const dimensioned<scalarMinMax> fMinMax
269  (
270  dimless,
271  scalarMinMax(Zero, scalar(1) - Foam::exp(-scalar(400)/scalar(50)))
272  );
273 
274  f.clip(fMinMax);
275 }
276 
277 
278 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
279 // Main program:
280 int main(int argc, char *argv[])
281 {
283  (
284  "Sets initial turbulence fields based on"
285  " various empirical equations"
286  );
288  argList::addOption("dict", "file", "Alternative setTurbulenceFieldsDict");
289 
290  #include "addRegionOption.H"
291 
292  #include "setRootCase.H"
293  #include "createTime.H"
294 
295  const word dictName("setTurbulenceFieldsDict");
297  Info<< "Reading " << dictIO.name() << nl << endl;
299 
300  #include "createNamedMesh.H"
301 
302  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
303 
305 
306 
307  // Dictionary input (M:p. 2)
308 
309  const scalar uRef = dict.getCheck<scalar>("uRef", scalarMinMax::ge(SMALL));
310  const dimensionedScalar uTau(dimVelocity, 0.05*uRef);
311  const scalar kappa =
312  dict.getCheckOrDefault<scalar>("kappa", 0.41, scalarMinMax::ge(SMALL));
313  const scalar Cmu =
314  dict.getCheckOrDefault<scalar>("Cmu", 0.09, scalarMinMax::ge(SMALL));
315  const scalar dPlusRef =
316  dict.getCheckOrDefault<scalar>("dPlusRef", 15, scalarMinMax::ge(SMALL));
317  const word fName = dict.getOrDefault<word>("f", "f");
318  const word UName = dict.getOrDefault<word>("U", "U");
319  const word epsilonName = dict.getOrDefault<word>("epsilon", "epsilon");
320  const word kName = dict.getOrDefault<word>("k", "k");
321  const word omegaName = dict.getOrDefault<word>("omega", "omega");
322  const word RName = dict.getOrDefault<word>("R", "R");
323  const bool initU = dict.getOrDefault<bool>("initialiseU", false);
324  const bool initEpsilon =
325  dict.getOrDefault<bool>("initialiseEpsilon", false);
326  const bool initK = dict.getOrDefault<bool>("initialiseK", false);
327  const bool initOmega = dict.getOrDefault<bool>("initialiseOmega", false);
328  const bool initR = dict.getOrDefault<bool>("initialiseR", false);
329  const bool writeF = dict.getOrDefault<bool>("writeF", false);
330 
331 
332  // Start initialising the operand fields
333 
334  // Read operand fields
335 
336  auto tU = tmp<volVectorField>::New
337  (
338  createIOobject(mesh, UName, IOobject::MUST_READ),
339  mesh
340  );
341 
342  // Infer the initial BCs from the velocity
343  const wordList bcTypes
344  (
345  tU.cref().boundaryField().size(),
346  fixedValueFvPatchScalarField::typeName
347  );
348 
349  tmp<volScalarField> tepsilon;
351  tmp<volScalarField> tomega;
353 
354  if (initEpsilon)
355  {
356  tepsilon = tmp<volScalarField>::New
357  (
358  createIOobject(mesh, epsilonName),
359  mesh,
361  bcTypes
362  );
363  }
364 
365  if (initK)
366  {
368  (
369  createIOobject(mesh, kName),
370  mesh,
372  bcTypes
373  );
374  }
375 
376  if (initOmega)
377  {
378  tomega = tmp<volScalarField>::New
379  (
380  createIOobject(mesh, omegaName),
381  mesh,
383  bcTypes
384  );
385  }
386 
387  if (initR)
388  {
390  (
391  createIOobject(mesh, RName),
392  mesh,
394  bcTypes
395  );
396  }
397 
398 
399  // Create elliptic blending factor field
400 
402  (
403  IOobject
404  (
405  fName,
406  runTime.timeName(),
407  mesh,
410  false // do not register
411  ),
412  mesh,
413  dimensionedScalar(dimless, scalar(1)),
415  );
416 
417  for (fvPatchScalarField& pfld : f.boundaryFieldRef())
418  {
419  if (isA<wallFvPatch>(pfld.patch()))
420  {
421  pfld == Zero;
422  }
423  }
424 
425 
426  // Create auxillary fields for the initialisation
427 
428  tmp<volScalarField> tnu = nu(mesh, tU.cref());
429  const volScalarField& nu = tnu.cref();
430 
431  // (M:p. 2)
432  tmp<volScalarField> tL = scalar(50)*nu/uTau;
433  const volScalarField& L = tL.cref();
434 
435  calcF(L, f);
436 
437  // (M:Eq. 8)
438  tmp<volScalarField> td = -tL*Foam::log(scalar(1) - f);
439  const volScalarField& d = td.cref();
440 
441  // (M:p. 2)
442  const volScalarField dPlus(d*uTau/nu);
443 
444  // (M:Eq. 11)
445  const volScalarField epsilon
446  (
447  pow4(uTau)/(kappa*nu*max(dPlus, dPlusRef))
448  );
449 
450  // (M:Eq. 13)
451  const volScalarField fk(Foam::exp(-dPlus/scalar(25)));
452 
453  // (M:Eq. 12)
454  const volScalarField k
455  (
456  (epsilon*sqr(td)*sqr(fk))/(2*nu)
457  + sqr(uTau)*sqr(scalar(1) - fk)/Foam::sqrt(Cmu)
458  );
459 
460 
461  // Initialise operand fields
462 
463  if (initU)
464  {
465  volVectorField& U = tU.ref();
466 
467  // Reichardt’s law (M:Eq. 10)
468  const scalar C = 7.8;
469  const scalar B1 = 11;
470  const scalar B2 = 3;
471  const volScalarField fRei
472  (
473  Foam::log(scalar(1) + kappa*dPlus)/kappa
474  + C*
475  (
476  scalar(1)
477  - Foam::exp(-dPlus/B1)
478  - dPlus/B1*Foam::exp(-dPlus/B2)
479  )
480  );
481 
482  // (M:Eq. 9)
483  const dimensionedScalar maxU(dimVelocity, SMALL);
484  U *= min(scalar(1), fRei*uTau/max(mag(U), maxU));
485  correctProcessorPatches<vector>(U);
486  }
487 
488  if (tepsilon.valid())
489  {
490  tepsilon.ref() = epsilon;
491  correctProcessorPatches<scalar>(tepsilon.ref());
492  }
493 
494  if (tk.valid())
495  {
496  tk.ref() = k;
497  correctProcessorPatches<scalar>(tk.ref());
498  }
499 
500  if (tomega.valid())
501  {
502  const dimensionedScalar k0(sqr(dimLength/dimTime), SMALL);
503  tomega.ref() = Cmu*epsilon/(k + k0);
504  correctProcessorPatches<scalar>(tomega.ref());
505  }
506 
507  if (tR.valid())
508  {
509  volSymmTensorField& R = tR.ref();
510 
511  // (M:Eq. 3)
512  const volSphericalTensorField Rdiag(k*twoThirdsI);
513  forAll(R, celli)
514  {
515  R[celli] = Rdiag[celli];
516  }
517  correctProcessorPatches<symmTensor>(R);
518  }
519 
520 
521  // Write operand fields
522 
523  Info<< endl;
524 
525  if (initU)
526  {
527  InfoField(tU->name());
528  tU->write();
529  }
530 
531  if (tepsilon.valid())
532  {
533  InfoField(tepsilon->name());
534  tepsilon->write();
535  }
536 
537  if (tk.valid())
538  {
539  InfoField(tk->name());
540  tk->write();
541  }
542 
543  if (tomega.valid())
544  {
545  InfoField(tomega->name());
546  tomega->write();
547  }
548 
549  if (tR.valid())
550  {
551  InfoField(tR->name());
552  tR->write();
553  }
554 
555  if (writeF)
556  {
557  InfoField(f.name());
558  f.write();
559  }
560 
561 
562  Info<< nl;
564 
565  Info<< "End\n" << endl;
566 
567  return 0;
568 }
569 
570 
571 // ************************************************************************* //
word dictName() const
The local dictionary name (final part of scoped name)
Definition: dictionaryI.H:53
Info<< "Reading thermophysical properties\"<< endl;autoPtr< psiReactionThermo > pThermo(psiReactionThermo::New(mesh))
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
bool valid() const noexcept
Identical to good(), or bool operator.
Definition: tmp.H:456
dictionary dict
static void noFunctionObjects(bool addWithOption=false)
Remove &#39;-noFunctionObjects&#39; option and ignore any occurrences.
Definition: argList.C:514
Graphite solid properties.
Definition: C.H:46
static void addNote(const string &note)
Add extra notes for the usage information.
Definition: argList.C:453
static const sphericalTensor twoThirdsI(2.0/3.0)
MinMax< scalar > scalarMinMax
A scalar min/max range.
Definition: MinMax.H:111
SolverPerformance< Type > solve(faMatrix< Type > &, const dictionary &)
Solve returning the solution statistics given convergence tolerance.
dimensionedScalar log(const dimensionedScalar &ds)
static autoPtr< fluidThermo > New(const fvMesh &, const word &phaseName=word::null)
Selector.
Definition: fluidThermo.C:68
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
const vector L(dict.get< vector >("L"))
const word UName(propsDict.getOrDefault< word >("U", "U"))
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:150
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition: tmpI.H:210
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
dimensionedSymmTensor sqr(const dimensionedVector &dv)
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:49
T getCheck(const word &keyword, const Predicate &pred, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T with additional checking FatalIOError if not found, or if the number of tokens is...
const word dictName("faMeshDefinition")
static autoPtr< IncompressibleTurbulenceModel > New(const volVectorField &U, const surfaceScalarField &phi, const TransportModel &transportModel, const word &propertiesName=turbulenceModel::propertiesName)
Return a reference to the selected turbulence model.
engineTime & runTime
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
const T & cref() const
Return const reference to the object or to the contents of a (non-null) managed pointer.
Definition: tmpI.H:196
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:639
Required Variables.
static unsigned int defaultPrecision() noexcept
Return the default precision.
Definition: IOstream.H:416
Generic GeometricField class.
Definition: areaFieldsFwd.H:50
label k
Boltzmann constant.
Generic dimensioned Type class.
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:361
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:50
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:413
scalar uTau
singlePhaseTransportModel laminarTransport(U, phi)
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:52
dimensionedScalar exp(const dimensionedScalar &ds)
A class for handling words, derived from Foam::string.
Definition: word.H:63
Fundamental fluid thermodynamic properties.
Definition: fluidThermo.H:48
static MinMax< scalar > ge(const scalar &minVal)
A semi-infinite range from minVal to the type max.
Definition: MinMaxI.H:24
T getCheckOrDefault(const word &keyword, const T &deflt, const Predicate &pred, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:203
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:376
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
static word timeName(const scalar t, const int precision=precision_)
Return time name of given scalar time formatted with the given precision.
Definition: Time.C:760
const word & constant() const noexcept
Return constant name.
Definition: TimePathsI.H:89
Info<< "Reading field U\"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho());volVectorField rhoU(IOobject("rhoU", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *U);volScalarField rhoE(IOobject("rhoE", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *(e+0.5 *magSqr(U)));surfaceScalarField pos(IOobject("pos", runTime.timeName(), mesh), mesh, dimensionedScalar("pos", dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar("neg", dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
Ostream & printExecutionTime(OSstream &os) const
Print the elapsed ExecutionTime (cpu-time), ClockTime.
Definition: TimeIO.C:612
labelList f(nPoints)
dimensioned< symmTensor > dimensionedSymmTensor
Dimensioned tensor obtained from generic dimensioned type.
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)
Return a reference to the selected turbulence model.
dimensionedScalar pow3(const dimensionedScalar &ds)
U
Definition: pEqn.H:72
#define R(A, B, C, D, E, F, K, M)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:41
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:79
scalar epsilon
dimensionedScalar pow4(const dimensionedScalar &ds)
Nothing to be read.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
messageStream Info
Information stream (stdout output on master, null elsewhere)
IOobject dictIO
This boundary condition enables processor communication across patches.
void clear() const noexcept
If object pointer points to valid object: delete object and set pointer to nullptr.
Definition: tmpI.H:263
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
const Type * isA(const U &obj)
Check if dynamic_cast to Type is possible.
Definition: typeInfo.H:87
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
A class for managing temporary objects.
Definition: HashPtrTable.H:50
virtual bool write(const bool valid=true) const
Write using setting from DB.
A simple single-phase transport model based on viscosityModel.
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:166
volScalarField & nu
Creates and initialises the face-flux field phi.
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:157
const dimensionSet dimVelocity
readOption
Enumeration defining read preferences.