turbulenceFields.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) 2013-2016 OpenFOAM Foundation
9  Copyright (C) 2015-2022 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 
29 #include "turbulenceFields.H"
32 #include "DESModelBase.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 namespace functionObjects
40 {
41  defineTypeNameAndDebug(turbulenceFields, 0);
42  addToRunTimeSelectionTable(functionObject, turbulenceFields, dictionary);
43 }
44 }
45 
46 const Foam::Enum
47 <
49 >
51 ({
52  { compressibleField::cfK, "k" },
53  { compressibleField::cfEpsilon, "epsilon" },
54  { compressibleField::cfOmega, "omega" },
55  { compressibleField::cfNuTilda, "nuTilda" },
56  { compressibleField::cfMut, "mut" },
57  { compressibleField::cfMuEff, "muEff" },
58  { compressibleField::cfAlphat, "alphat" },
59  { compressibleField::cfAlphaEff, "alphaEff" },
60  { compressibleField::cfR, "R" },
61  { compressibleField::cfDevRhoReff, "devRhoReff" },
62  { compressibleField::cfL, "L" },
63  { compressibleField::cfI, "I" },
64  { compressibleField::cfLESRegion, "LESRegion" },
65  { compressibleField::cffd, "fd" },
66 });
67 
68 
69 const Foam::Enum
70 <
72 >
74 ({
75  { incompressibleField::ifK, "k" },
76  { incompressibleField::ifEpsilon, "epsilon" },
77  { incompressibleField::ifOmega, "omega" },
78  { incompressibleField::ifNuTilda, "nuTilda" },
79  { incompressibleField::ifNut, "nut" },
80  { incompressibleField::ifNuEff, "nuEff" },
81  { incompressibleField::ifR, "R" },
82  { incompressibleField::ifDevReff, "devReff" },
83  { incompressibleField::ifL, "L" },
84  { incompressibleField::ifI, "I" },
85  { incompressibleField::ifLESRegion, "LESRegion" },
86  { incompressibleField::iffd, "fd" },
87 });
88 
89 
91 (
93 );
94 
95 
96 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
97 
99 {
100  for (const word& f : fieldSet_)
101  {
102  const word localName(IOobject::scopedName(prefix_, f));
103 
104  if (obr_.found(localName))
105  {
107  << "Cannot store turbulence field " << localName
108  << " since an object with that name already exists"
109  << nl << endl;
110 
111  fieldSet_.unset(f);
112  }
113  }
114 
115  initialised_ = true;
116 }
117 
118 
120 {
121  if (obr_.foundObject<compressible::turbulenceModel>(modelName_))
122  {
123  return true;
124  }
125  else if (obr_.foundObject<incompressible::turbulenceModel>(modelName_))
126  {
127  return false;
128  }
129 
131  << "Turbulence model not found in database, deactivating"
132  << exit(FatalError);
133 
134  return false;
135 }
136 
137 
138 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
139 
141 (
142  const word& name,
143  const Time& runTime,
144  const dictionary& dict
145 )
146 :
148  initialised_(false),
149  prefix_(dict.getOrDefault<word>("prefix", "turbulenceProperties")),
150  fieldSet_()
151 {
152  read(dict);
153 }
154 
155 
156 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
157 
159 {
161  {
162  dict.readIfPresent("prefix", prefix_);
163 
164  if (dict.found("field"))
165  {
166  fieldSet_.insert(dict.get<word>("field"));
167  }
168  else
169  {
170  fieldSet_.insert(dict.get<wordList>("fields"));
171  }
172 
173  Info<< type() << " " << name() << ": ";
174  if (fieldSet_.size())
175  {
176  Info<< "storing fields:" << nl;
177  for (const word& f : fieldSet_)
178  {
179  Info<< " " << IOobject::scopedName(prefix_, f) << nl;
180  }
181  Info<< endl;
182  }
183  else
184  {
185  Info<< "no fields requested to be stored" << nl << endl;
186  }
187 
188  initialised_ = false;
189 
190  return true;
191  }
192 
193  return false;
194 }
195 
196 
198 {
199  if (!initialised_)
200  {
201  initialise();
202  }
203 
204  const bool comp = compressible();
205 
206  if (comp)
207  {
208  const auto& model =
209  obr_.lookupObject<compressible::turbulenceModel>(modelName_);
210 
211  for (const word& f : fieldSet_)
212  {
213  switch (compressibleFieldNames_[f])
214  {
215  case cfK:
216  {
217  processField<scalar>(f, model.k());
218  break;
219  }
220  case cfEpsilon:
221  {
222  processField<scalar>(f, model.epsilon());
223  break;
224  }
225  case cfOmega:
226  {
227  processField<scalar>(f, model.omega());
228  break;
229  }
230  case cfNuTilda:
231  {
232  processField<scalar>(f, nuTilda(model));
233  break;
234  }
235  case cfMut:
236  {
237  processField<scalar>(f, model.mut());
238  break;
239  }
240  case cfMuEff:
241  {
242  processField<scalar>(f, model.muEff());
243  break;
244  }
245  case cfAlphat:
246  {
247  processField<scalar>(f, model.alphat());
248  break;
249  }
250  case cfAlphaEff:
251  {
252  processField<scalar>(f, model.alphaEff());
253  break;
254  }
255  case cfR:
256  {
257  processField<symmTensor>(f, model.R());
258  break;
259  }
260  case cfDevRhoReff:
261  {
262  processField<symmTensor>(f, model.devRhoReff());
263  break;
264  }
265  case cfL:
266  {
267  processField<scalar>(f, L(model));
268  break;
269  }
270  case cfI:
271  {
272  processField<scalar>(f, I(model));
273  break;
274  }
275  case cfLESRegion:
276  {
277  auto* DESPtr = mesh_.cfindObject<DESModelBase>(modelName_);
278  if (!DESPtr)
279  {
281  << "Turbulence model is not a DES model - "
282  << "skipping request for LESRegion" << endl;
283 
284  break;
285  }
286 
287  processField<scalar>(f, DESPtr->LESRegion());
288  break;
289  }
290  case cffd:
291  {
292  auto* DESPtr = mesh_.cfindObject<DESModelBase>(modelName_);
293  if (!DESPtr)
294  {
296  << "Turbulence model is not a DES model - "
297  << "skipping request for fd" << endl;
298 
299  break;
300  }
301 
302  processField<scalar>(f, DESPtr->fd());
303  break;
304  }
305  default:
306  {
308  << "Invalid field selection" << abort(FatalError);
309  }
310  }
311  }
312  }
313  else
314  {
315  const auto& model =
316  obr_.lookupObject<incompressible::turbulenceModel>(modelName_);
317 
318  for (const word& f : fieldSet_)
319  {
320  switch (incompressibleFieldNames_[f])
321  {
322  case ifK:
323  {
324  processField<scalar>(f, model.k());
325  break;
326  }
327  case ifEpsilon:
328  {
329  processField<scalar>(f, model.epsilon());
330  break;
331  }
332  case ifOmega:
333  {
334  processField<scalar>(f, model.omega());
335  break;
336  }
337  case ifNuTilda:
338  {
339  processField<scalar>(f, nuTilda(model));
340  break;
341  }
342  case ifNut:
343  {
344  processField<scalar>(f, model.nut());
345  break;
346  }
347  case ifNuEff:
348  {
349  processField<scalar>(f, model.nuEff());
350  break;
351  }
352  case ifR:
353  {
354  processField<symmTensor>(f, model.R());
355  break;
356  }
357  case ifDevReff:
358  {
359  processField<symmTensor>(f, model.devReff());
360  break;
361  }
362  case ifL:
363  {
364  processField<scalar>(f, L(model));
365  break;
366  }
367  case ifI:
368  {
369  processField<scalar>(f, I(model));
370  break;
371  }
372  case ifLESRegion:
373  {
374  auto* DESPtr = mesh_.cfindObject<DESModelBase>(modelName_);
375  if (!DESPtr)
376  {
378  << "Turbulence model is not a DES model - "
379  << "skipping request for LESRegion" << endl;
380 
381  break;
382  }
383 
384  processField<scalar>(f, DESPtr->LESRegion());
385  break;
386  }
387  case iffd:
388  {
389  auto* DESPtr = mesh_.cfindObject<DESModelBase>(modelName_);
390  if (!DESPtr)
391  {
393  << "Turbulence model is not a DES model - "
394  << "skipping request for fd" << endl;
395 
396  break;
397  }
398 
399  processField<scalar>(f, DESPtr->fd());
400  break;
401  }
402  default:
403  {
405  << "Invalid field selection" << abort(FatalError);
406  }
407  }
408  }
409  }
410 
411  return true;
412 }
413 
414 
416 {
417  for (const word& f : fieldSet_)
418  {
419  const word localName(IOobject::scopedName(prefix_, f));
420 
421  writeObject(localName);
422  }
423  Info<< endl;
424 
425  return true;
426 }
427 
428 
429 // ************************************************************************* //
dictionary dict
defineTypeNameAndDebug(ObukhovLength, 0)
static bool initialised_(false)
void initialise()
Unset duplicate fields already registered by other function objects.
bool initialised_
Flag to track initialisation.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
compressibleField
Options for the turbulence fields (compressible)
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
const vector L(dict.get< vector >("L"))
turbulenceFields(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
engineTime & runTime
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
static const word modelName_
Name of the turbulence properties dictionary.
IncompressibleTurbulenceModel< transportModel > turbulenceModel
word prefix_
Name of output-field prefix.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
Macros for easy insertion into run-time selection tables.
wordHashSet fieldSet_
Fields to load.
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:127
bool compressible()
Return true if compressible turbulence model is identified.
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition: IOobjectI.H:40
bool unset(const Key &key)
Unset the specified key - same as erase.
Definition: HashSet.H:250
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
static const Identity< scalar > I
Definition: Identity.H:100
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition: word.H:63
virtual bool read(const dictionary &)
Read the controls.
errorManip< error > abort(error &err)
Definition: errorManip.H:139
static const Enum< incompressibleField > incompressibleFieldNames_
Names for incompressibleField turbulence fields.
static const Enum< compressibleField > compressibleFieldNames_
Names for compressibleField turbulence fields.
bool compressible
Definition: pEqn.H:2
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
labelList f(nPoints)
#define WarningInFunction
Report a warning using Foam::Warning.
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: error.H:64
const objectRegistry & obr_
Reference to the region objectRegistry.
virtual bool execute()
Calculate turbulence fields.
bool found(const word &name, bool recursive=false) const
Same as contains()
messageStream Info
Information stream (stdout output on master, null elsewhere)
virtual bool read(const dictionary &dict)
Read optional controls.
incompressibleField
Options for the turbulence fields (incompressible)
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
Namespace for OpenFOAM.