turbulenceFields.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2013-2018 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "turbulenceFields.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace functionObjects
36 {
37  defineTypeNameAndDebug(turbulenceFields, 0);
38 
40  (
41  functionObject,
42  turbulenceFields,
43  dictionary
44  );
45 }
46 }
47 
48 template<>
49 const char* Foam::NamedEnum
50 <
52  9
53 >::names[] =
54 {
55  "k",
56  "epsilon",
57  "omega",
58  "mut",
59  "muEff",
60  "alphat",
61  "alphaEff",
62  "R",
63  "devRhoReff"
64 };
65 
66 const Foam::NamedEnum
67 <
69  9
71 
72 template<>
73 const char* Foam::NamedEnum
74 <
76  7
77 >::names[] =
78 {
79  "k",
80  "epsilon",
81  "omega",
82  "nut",
83  "nuEff",
84  "R",
85  "devReff"
86 };
87 
88 const Foam::NamedEnum
89 <
91  7
93 
94 
95 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
96 
98 (
99  const word& name,
100  const Time& runTime,
101  const dictionary& dict
102 )
103 :
104  fvMeshFunctionObject(name, runTime, dict),
105  fieldSet_()
106 {
107  read(dict);
108 }
109 
110 
111 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
112 
114 {}
115 
116 
117 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
118 
120 {
122 }
123 
124 
126 {
127  if (dict.found("field"))
128  {
129  fieldSet_.insert(word(dict.lookup("field")));
130  }
131  else
132  {
133  fieldSet_.insert(wordList(dict.lookup("fields")));
134  }
135 
136  if (dict.lookupOrDefault<Switch>("prefix", false))
137  {
138  prefix_ = modelName() + ':';
139  }
140  else
141  {
142  prefix_ = word::null;
143  }
144 
145  Info<< type() << " " << name() << ": ";
146  if (fieldSet_.size())
147  {
148  Info<< "storing fields:" << nl;
149  forAllConstIter(wordHashSet, fieldSet_, iter)
150  {
151  Info<< " " << prefix_ + iter.key() << nl;
152  }
153  Info<< endl;
154  }
155  else
156  {
157  Info<< "no fields requested to be stored" << nl << endl;
158  }
159 
160  return true;
161 }
162 
163 
165 {
166  if (obr_.foundObject<compressible::turbulenceModel>(modelName()))
167  {
168  const compressible::turbulenceModel& model =
169  obr_.lookupObject<compressible::turbulenceModel>(modelName());
170 
171  forAllConstIter(wordHashSet, fieldSet_, iter)
172  {
173  const word& f = iter.key();
174  switch (compressibleFieldNames_[f])
175  {
177  {
178  processField<scalar>(f, model.k());
179  break;
180  }
182  {
183  processField<scalar>(f, model.epsilon());
184  break;
185  }
186  case compressibleField::omega:
187  {
188  processField<scalar>(f, omega(model));
189  break;
190  }
191  case compressibleField::mut:
192  {
193  processField<scalar>(f, model.mut());
194  break;
195  }
196  case compressibleField::muEff:
197  {
198  processField<scalar>(f, model.muEff());
199  break;
200  }
201  case compressibleField::alphat:
202  {
203  processField<scalar>(f, model.alphat());
204  break;
205  }
206  case compressibleField::alphaEff:
207  {
208  processField<scalar>(f, model.alphaEff());
209  break;
210  }
212  {
213  processField<symmTensor>(f, model.R());
214  break;
215  }
216  case compressibleField::devRhoReff:
217  {
218  processField<symmTensor>(f, model.devRhoReff());
219  break;
220  }
221  default:
222  {
224  << "Invalid field selection" << exit(FatalError);
225  }
226  }
227  }
228  }
229  else if (obr_.foundObject<compressibleTurbulenceModel>(modelName()))
230  {
231  const compressibleTurbulenceModel& model =
232  obr_.lookupObject<compressibleTurbulenceModel>(modelName());
233 
234  forAllConstIter(wordHashSet, fieldSet_, iter)
235  {
236  const word& f = iter.key();
237  switch (compressibleFieldNames_[f])
238  {
240  {
241  processField<scalar>(f, model.k());
242  break;
243  }
245  {
246  processField<scalar>(f, model.epsilon());
247  break;
248  }
249  case compressibleField::omega:
250  {
251  processField<scalar>(f, omega(model));
252  break;
253  }
254  case compressibleField::mut:
255  {
256  processField<scalar>(f, model.mut());
257  break;
258  }
259  case compressibleField::muEff:
260  {
261  processField<scalar>(f, model.muEff());
262  break;
263  }
265  {
266  processField<symmTensor>(f, model.R());
267  break;
268  }
269  case compressibleField::devRhoReff:
270  {
271  processField<symmTensor>(f, model.devRhoReff());
272  break;
273  }
274  default:
275  {
277  << "Invalid field selection" << exit(FatalError);
278  }
279  }
280  }
281  }
282  else if (obr_.foundObject<incompressible::turbulenceModel>(modelName()))
283  {
284  const incompressible::turbulenceModel& model =
285  obr_.lookupObject<incompressible::turbulenceModel>(modelName());
286 
287  forAllConstIter(wordHashSet, fieldSet_, iter)
288  {
289  const word& f = iter.key();
290  switch (incompressibleFieldNames_[f])
291  {
293  {
294  processField<scalar>(f, model.k());
295  break;
296  }
298  {
299  processField<scalar>(f, model.epsilon());
300  break;
301  }
302  case incompressibleField::omega:
303  {
304  processField<scalar>(f, omega(model));
305  break;
306  }
308  {
309  processField<scalar>(f, model.nut());
310  break;
311  }
312  case incompressibleField::nuEff:
313  {
314  processField<scalar>(f, model.nuEff());
315  break;
316  }
318  {
319  processField<symmTensor>(f, model.R());
320  break;
321  }
322  case incompressibleField::devReff:
323  {
324  processField<symmTensor>(f, model.devReff());
325  break;
326  }
327  default:
328  {
330  << "Invalid field selection" << exit(FatalError);
331  }
332  }
333  }
334  }
335  else
336  {
338  << "Turbulence model not found in database, deactivating"
339  << exit(FatalError);
340  }
341 
342  return true;
343 }
344 
345 
347 {
348  forAllConstIter(wordHashSet, fieldSet_, iter)
349  {
350  const word fieldName = prefix_ + iter.key();
351  writeObject(fieldName);
352  }
353 
354  return true;
355 }
356 
357 
358 // ************************************************************************* //
A HashTable with keys but without contents.
Definition: HashSet.H:59
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:438
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
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:319
addToRunTimeSelectionTable(functionObject, Qdot, dictionary)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, y/n, t/f, or none/any.
Definition: Switch.H:60
virtual tmp< volScalarField > muEff() const =0
Return the effective dynamic viscosity.
virtual tmp< volSymmTensorField > R() const =0
Return the Reynolds stress tensor.
virtual tmp< volScalarField > mut() const =0
Return the turbulence dynamic viscosity.
static const NamedEnum< incompressibleField, 7 > incompressibleFieldNames_
label k
Boltzmann constant.
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:51
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
Macros for easy insertion into run-time selection tables.
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
virtual tmp< volSymmTensorField > devRhoReff() const =0
Return the effective stress tensor including the laminar stress.
virtual tmp< volScalarField > alphaEff() const
Effective thermal turbulent diffusivity of mixture [kg/m/s].
Templated abstract base class for single-phase incompressible turbulence models.
bool read(const char *, int32_t &)
Definition: int32IO.C:85
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from string.
Definition: word.H:59
static const word null
An empty word.
Definition: word.H:77
virtual bool read(const dictionary &)
Read the controls.
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
virtual tmp< volSymmTensorField > devReff() const
Return the effective stress tensor.
static const char nl
Definition: Ostream.H:265
Abstract base class for turbulence models (RAS, LES and laminar).
virtual tmp< volScalarField > nuEff() const =0
Return the effective viscosity.
labelList f(nPoints)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
defineTypeNameAndDebug(Qdot, 0)
List< word > wordList
A List of words.
Definition: fileName.H:54
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
#define R(A, B, C, D, E, F, K, M)
virtual tmp< volScalarField > epsilon() const =0
Return the turbulence kinetic energy dissipation rate.
scalar epsilon
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
virtual bool execute()
Calculate turbulence fields.
messageStream Info
virtual tmp< volScalarField > alphat() const
Turbulent thermal diffusivity for enthalpy [kg/m/s].
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
virtual tmp< volScalarField > nut() const =0
Return the turbulence viscosity.
static const NamedEnum< compressibleField, 9 > compressibleFieldNames_
scalar nut
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:583