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-2025 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 {
38 
40  (
44  );
45 }
46 }
47 
48 const Foam::NamedEnum
49 <
51  7
53 {
54  "k",
55  "epsilon",
56  "omega",
57  "nut",
58  "nuEff",
59  "kappaEff",
60  "R"
61 };
62 
63 const Foam::NamedEnum
64 <
66  6
68 {
69  "k",
70  "epsilon",
71  "omega",
72  "nut",
73  "nuEff",
74  "R"
75 };
76 
77 
78 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
79 
81 (
82  const word& name,
83  const Time& runTime,
84  const dictionary& dict
85 )
86 :
87  fvMeshFunctionObject(name, runTime, dict),
88  fieldSet_(),
89  phaseName_(dict.lookupOrDefault<word>("phase", word::null))
90 {
91  read(dict);
92 }
93 
94 
95 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
96 
98 {}
99 
100 
101 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
102 
104 {
105  if (dict.found("field"))
106  {
107  fieldSet_.insert(word(dict.lookup("field")));
108  }
109  else
110  {
111  fieldSet_.insert(wordList(dict.lookup("fields")));
112  }
113 
114  if (dict.lookupOrDefault<Switch>("prefix", false))
115  {
116  prefix_ = momentumTransportModel::typeName + ':';
117  }
118  else
119  {
120  prefix_ = word::null;
121  }
122 
123  Info<< type() << " " << name() << ": ";
124  if (fieldSet_.size())
125  {
126  Info<< "storing fields:" << nl;
127  forAllConstIter(wordHashSet, fieldSet_, iter)
128  {
129  Info<< " "
130  << IOobject::groupName(prefix_ + iter.key(), phaseName_) << nl;
131  }
132  Info<< endl;
133  }
134  else
135  {
136  Info<< "no fields requested to be stored" << nl << endl;
137  }
138 
139  return true;
140 }
141 
142 
144 {
145  if (obr_.foundType<fluidThermophysicalTransportModel>(phaseName_))
146  {
148  obr_.lookupType<fluidThermophysicalTransportModel>(phaseName_);
149 
151  ttm.momentumTransport();
152 
153  forAllConstIter(wordHashSet, fieldSet_, iter)
154  {
155  const word& f = iter.key();
156  switch (compressibleFieldNames_[f])
157  {
159  {
160  processField<scalar>(f, model.k());
161  break;
162  }
164  {
165  processField<scalar>(f, model.epsilon());
166  break;
167  }
169  {
170  processField<scalar>(f, model.omega());
171  break;
172  }
174  {
175  processField<scalar>(f, model.nut());
176  break;
177  }
179  {
180  processField<scalar>(f, model.nuEff());
181  break;
182  }
183  case compressibleField::kappaEff:
184  {
185  processField<scalar>(f, ttm.kappaEff());
186  break;
187  }
189  {
190  processField<symmTensor>(f, model.R());
191  break;
192  }
193  default:
194  {
196  << "Invalid field selection" << exit(FatalError);
197  }
198  }
199  }
200  }
201  else if (obr_.foundType<compressibleMomentumTransportModel>(phaseName_))
202  {
204  obr_.lookupType<compressibleMomentumTransportModel>(phaseName_);
205 
206  forAllConstIter(wordHashSet, fieldSet_, iter)
207  {
208  const word& f = iter.key();
209  switch (compressibleFieldNames_[f])
210  {
212  {
213  processField<scalar>(f, model.k());
214  break;
215  }
217  {
218  processField<scalar>(f, model.epsilon());
219  break;
220  }
222  {
223  processField<scalar>(f, model.omega());
224  break;
225  }
227  {
228  processField<scalar>(f, model.nut());
229  break;
230  }
232  {
233  processField<scalar>(f, model.nuEff());
234  break;
235  }
237  {
238  processField<symmTensor>(f, model.R());
239  break;
240  }
241  default:
242  {
244  << "Invalid field selection" << exit(FatalError);
245  }
246  }
247  }
248  }
249  else if
250  (
251  obr_.foundType<incompressible::momentumTransportModel>(phaseName_)
252  )
253  {
255  obr_.lookupType<incompressible::momentumTransportModel>(phaseName_);
256 
257  forAllConstIter(wordHashSet, fieldSet_, iter)
258  {
259  const word& f = iter.key();
260  switch (incompressibleFieldNames_[f])
261  {
263  {
264  processField<scalar>(f, model.k());
265  break;
266  }
268  {
269  processField<scalar>(f, model.epsilon());
270  break;
271  }
273  {
274  processField<scalar>(f, model.omega());
275  break;
276  }
278  {
279  processField<scalar>(f, model.nut());
280  break;
281  }
283  {
284  processField<scalar>(f, model.nuEff());
285  break;
286  }
288  {
289  processField<symmTensor>(f, model.R());
290  break;
291  }
292  default:
293  {
295  << "Invalid field selection" << exit(FatalError);
296  }
297  }
298  }
299  }
300  else
301  {
303  << "Turbulence model not found in database, deactivating"
304  << exit(FatalError);
305  }
306 
307  return true;
308 }
309 
310 
312 {
313  forAllConstIter(wordHashSet, fieldSet_, iter)
314  {
315  const word fieldName
316  (
317  IOobject::groupName(prefix_ + iter.key(), phaseName_)
318  );
319  writeObject(fieldName);
320  }
321 
322  return true;
323 }
324 
325 
326 // ************************************************************************* //
label k
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:476
Macros for easy insertion into run-time selection tables.
A HashTable with keys but without contents.
Definition: HashSet.H:62
static word groupName(Name name, const word &group)
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:55
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:61
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:76
Base class for single-phase compressible turbulence models.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Abstract base class for fluid thermophysical transport models RAS, LES and laminar.
const compressibleMomentumTransportModel & momentumTransport() const
Access function to momentum transport model.
Abstract base-class for Time/database functionObjects.
Specialisation of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
Stores derived turbulence fields on the mesh database for further manipulation.
turbulenceFields(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
static const NamedEnum< incompressibleField, 6 > incompressibleFieldNames_
static const NamedEnum< compressibleField, 7 > compressibleFieldNames_
virtual bool execute()
Calculate turbulence fields.
virtual bool write()
Write the turbulence fields.
virtual bool read(const dictionary &)
Read the controls.
Base class for single-phase incompressible turbulence models.
virtual tmp< volScalarField > nut() const =0
Return the turbulence viscosity.
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
virtual tmp< volSymmTensorField > R() const =0
Return the Reynolds stress tensor [m^2/s^2].
virtual tmp< volScalarField > omega() const =0
Return the turbulence specific dissipation rate.
virtual tmp< volScalarField > epsilon() const =0
Return the turbulence kinetic energy dissipation rate.
virtual tmp< volScalarField > nuEff() const =0
Return the effective viscosity.
virtual tmp< volScalarField > kappaEff() const =0
Effective thermal turbulent conductivity.
A class for handling words, derived from string.
Definition: word.H:62
static const word null
An empty word.
Definition: word.H:77
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
const scalar omega
const scalar nut
const scalar nuEff
const scalar epsilon
const volSymmTensorField R(IOobject("R", runTime.name(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), turbulence->R())
defineTypeNameAndDebug(adjustTimeStepToCombustion, 0)
addToRunTimeSelectionTable(functionObject, adjustTimeStepToCombustion, dictionary)
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
List< word > wordList
A List of words.
Definition: fileName.H:54
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
messageStream Info
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
error FatalError
static const char nl
Definition: Ostream.H:267
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
labelList f(nPoints)
dictionary dict