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-2020 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"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 namespace functionObjects
39 {
40  defineTypeNameAndDebug(turbulenceFields, 0);
41  addToRunTimeSelectionTable(functionObject, turbulenceFields, dictionary);
42 }
43 }
44 
45 const Foam::Enum
46 <
48 >
50 ({
51  { compressibleField::cfK, "k" },
52  { compressibleField::cfEpsilon, "epsilon" },
53  { compressibleField::cfOmega, "omega" },
54  { compressibleField::cfNuTilda, "nuTilda" },
55  { compressibleField::cfMut, "mut" },
56  { compressibleField::cfMuEff, "muEff" },
57  { compressibleField::cfAlphat, "alphat" },
58  { compressibleField::cfAlphaEff, "alphaEff" },
59  { compressibleField::cfR, "R" },
60  { compressibleField::cfDevRhoReff, "devRhoReff" },
61  { compressibleField::cfL, "L" },
62  { compressibleField::cfI, "I" },
63 });
64 
65 
66 const Foam::Enum
67 <
69 >
71 ({
72  { incompressibleField::ifK, "k" },
73  { incompressibleField::ifEpsilon, "epsilon" },
74  { incompressibleField::ifOmega, "omega" },
75  { incompressibleField::ifNuTilda, "nuTilda" },
76  { incompressibleField::ifNut, "nut" },
77  { incompressibleField::ifNuEff, "nuEff" },
78  { incompressibleField::ifR, "R" },
79  { incompressibleField::ifDevReff, "devReff" },
80  { incompressibleField::ifL, "L" },
81  { incompressibleField::ifI, "I" },
82 });
83 
84 
86 (
88 );
89 
90 
91 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
92 
94 {
96  {
97  return true;
98  }
100  {
101  return false;
102  }
103 
105  << "Turbulence model not found in database, deactivating"
106  << exit(FatalError);
107 
108  return false;
109 }
110 
111 
112 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
113 
115 (
116  const word& name,
117  const Time& runTime,
118  const dictionary& dict
119 )
120 :
122  fieldSet_()
123 {
124  read(dict);
125 }
126 
127 
128 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
129 
131 {
133 
134  if (dict.found("field"))
135  {
136  fieldSet_.insert(dict.get<word>("field"));
137  }
138  else
139  {
140  fieldSet_.insert(dict.get<wordList>("fields"));
141  }
142 
143  Info<< type() << " " << name() << ": ";
144  if (fieldSet_.size())
145  {
146  Info<< "storing fields:" << nl;
147  for (const word& f : fieldSet_)
148  {
149  Info<< " " << modelName_ << ':' << f << nl;
150  }
151  Info<< endl;
152  }
153  else
154  {
155  Info<< "no fields requested to be stored" << nl << endl;
156  }
157 
158  return true;
159 }
160 
161 
163 {
164  bool comp = compressible();
165 
166  if (comp)
167  {
168  const compressible::turbulenceModel& model =
169  obr_.lookupObject<compressible::turbulenceModel>(modelName_);
170 
171  for (const word& f : fieldSet_)
172  {
173  switch (compressibleFieldNames_[f])
174  {
175  case cfK:
176  {
177  processField<scalar>(f, model.k());
178  break;
179  }
180  case cfEpsilon:
181  {
182  processField<scalar>(f, model.epsilon());
183  break;
184  }
185  case cfOmega:
186  {
187  processField<scalar>(f, omega(model));
188  break;
189  }
190  case cfNuTilda:
191  {
192  processField<scalar>(f, nuTilda(model));
193  break;
194  }
195  case cfMut:
196  {
197  processField<scalar>(f, model.mut());
198  break;
199  }
200  case cfMuEff:
201  {
202  processField<scalar>(f, model.muEff());
203  break;
204  }
205  case cfAlphat:
206  {
207  processField<scalar>(f, model.alphat());
208  break;
209  }
210  case cfAlphaEff:
211  {
212  processField<scalar>(f, model.alphaEff());
213  break;
214  }
215  case cfR:
216  {
217  processField<symmTensor>(f, model.R());
218  break;
219  }
220  case cfDevRhoReff:
221  {
222  processField<symmTensor>(f, model.devRhoReff());
223  break;
224  }
225  case cfL:
226  {
227  processField<scalar>(f, L(model));
228  break;
229  }
230  case cfI:
231  {
232  processField<scalar>(f, I(model));
233  break;
234  }
235  default:
236  {
238  << "Invalid field selection" << abort(FatalError);
239  }
240  }
241  }
242  }
243  else
244  {
245  const incompressible::turbulenceModel& model =
246  obr_.lookupObject<incompressible::turbulenceModel>(modelName_);
247 
248  for (const word& f : fieldSet_)
249  {
250  switch (incompressibleFieldNames_[f])
251  {
252  case ifK:
253  {
254  processField<scalar>(f, model.k());
255  break;
256  }
257  case ifEpsilon:
258  {
259  processField<scalar>(f, model.epsilon());
260  break;
261  }
262  case ifOmega:
263  {
264  processField<scalar>(f, omega(model));
265  break;
266  }
267  case ifNuTilda:
268  {
269  processField<scalar>(f, nuTilda(model));
270  break;
271  }
272  case ifNut:
273  {
274  processField<scalar>(f, model.nut());
275  break;
276  }
277  case ifNuEff:
278  {
279  processField<scalar>(f, model.nuEff());
280  break;
281  }
282  case ifR:
283  {
284  processField<symmTensor>(f, model.R());
285  break;
286  }
287  case ifDevReff:
288  {
289  processField<symmTensor>(f, model.devReff());
290  break;
291  }
292  case ifL:
293  {
294  processField<scalar>(f, L(model));
295  break;
296  }
297  case ifI:
298  {
299  processField<scalar>(f, I(model));
300  break;
301  }
302  default:
303  {
305  << "Invalid field selection" << abort(FatalError);
306  }
307  }
308  }
309  }
310 
311  return true;
312 }
313 
314 
316 {
317  for (const word& f : fieldSet_)
318  {
319  const word fieldName = modelName_ + ':' + f;
320  writeObject(fieldName);
321  }
322 
323  return true;
324 }
325 
326 
327 // ************************************************************************* //
Foam::turbulenceModel::R
virtual tmp< volSymmTensorField > R() const =0
Return the Reynolds stress tensor.
runTime
engineTime & runTime
Definition: createEngineTime.H:13
L
const vector L(dict.get< vector >("L"))
Foam::functionObjects::regionFunctionObject::obr_
const objectRegistry & obr_
Reference to the region objectRegistry.
Definition: regionFunctionObject.H:103
Foam::Enum
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: IOstreamOption.H:57
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:73
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::functionObjects::turbulenceFields::incompressibleFieldNames_
static const Enum< incompressibleField > incompressibleFieldNames_
Names for incompressibleField turbulence fields.
Definition: turbulenceFields.H:230
Foam::turbulenceModel::nut
virtual tmp< volScalarField > nut() const =0
Return the turbulence viscosity.
Foam::functionObjects::turbulenceFields::modelName_
static const word modelName_
Turbulence closure model name.
Definition: turbulenceFields.H:233
Foam::read
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
Foam::functionObjects::turbulenceFields::compressibleField
compressibleField
Options for the turbulence fields (compressible)
Definition: turbulenceFields.H:195
turbulentTransportModel.H
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:100
Foam::functionObjects::turbulenceFields::turbulenceFields
turbulenceFields(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
Definition: turbulenceFields.C:115
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::functionObjects::fvMeshFunctionObject
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
Definition: fvMeshFunctionObject.H:64
Foam::objectRegistry::foundObject
bool foundObject(const word &name, const bool recursive=false) const
Is the named Type found?
Definition: objectRegistryTemplates.C:379
Foam::IncompressibleTurbulenceModel::devReff
virtual tmp< volSymmTensorField > devReff() const
Return the effective stress tensor.
Definition: IncompressibleTurbulenceModel.C:104
Foam::functionObjects::turbulenceFields::write
virtual bool write()
Do nothing.
Definition: turbulenceFields.C:315
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
compressible
bool compressible
Definition: pEqn.H:2
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::functionObjects::regionFunctionObject::read
virtual bool read(const dictionary &dict)
Read optional controls.
Definition: regionFunctionObject.C:173
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::turbulenceModel::k
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::functionObjects::turbulenceFields::read
virtual bool read(const dictionary &)
Read the controls.
Definition: turbulenceFields.C:130
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::functionObjects::turbulenceFields::compressibleFieldNames_
static const Enum< compressibleField > compressibleFieldNames_
Names for compressibleField turbulence fields.
Definition: turbulenceFields.H:212
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
Foam::nl
constexpr char nl
Definition: Ostream.H:385
Foam::functionObjects::turbulenceFields::incompressibleField
incompressibleField
Options for the turbulence fields (incompressible)
Definition: turbulenceFields.H:215
Foam::ThermalDiffusivity::alphaEff
virtual tmp< volScalarField > alphaEff() const
Return the effective turbulent thermal diffusivity for enthalpy.
Definition: ThermalDiffusivity.H:160
f
labelList f(nPoints)
Foam::functionObjects::addToRunTimeSelectionTable
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
Foam::ThermalDiffusivity
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
Definition: phaseCompressibleTurbulenceModelFwd.H:47
Foam::ThermalDiffusivity::alphat
virtual tmp< volScalarField > alphat() const
Return the turbulent thermal diffusivity for enthalpy [kg/m/s].
Definition: ThermalDiffusivity.C:121
Foam::List< word >
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::functionObjects::defineTypeNameAndDebug
defineTypeNameAndDebug(ObukhovLength, 0)
turbulenceFields.H
Foam::IncompressibleTurbulenceModel
Templated abstract base class for single-phase incompressible turbulence models.
Definition: IncompressibleTurbulenceModel.H:55
Foam::functionObjects::turbulenceFields::compressible
bool compressible()
Return true if compressible turbulence model is identified.
Definition: turbulenceFields.C:93
Foam::functionObjects::turbulenceFields::execute
virtual bool execute()
Calculate turbulence fields.
Definition: turbulenceFields.C:162
Foam::turbulenceModel::nuEff
virtual tmp< volScalarField > nuEff() const =0
Return the effective viscosity.
Foam::turbulenceModel::epsilon
virtual tmp< volScalarField > epsilon() const =0
Return the turbulence kinetic energy dissipation rate.
turbulentFluidThermoModel.H
Foam::I
static const Identity< scalar > I
Definition: Identity.H:95