turbulenceFields.H
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-2021 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 Class
28  Foam::functionObjects::turbulenceFields
29 
30 Group
31  grpFieldFunctionObjects
32 
33 Description
34  Computes various turbulence-related quantities that are not typically
35  output during calculations, and stores/writes them on the mesh database
36  for further manipulation.
37 
38  Fields are stored as copies of the original with a user-defined prefix
39  e.g. a prefix of \c turbulenceModel yields the following for field \c R:
40 
41  \verbatim
42  turbulenceModel:R
43  \endverbatim
44 
45  Operands:
46  \table
47  Operand | Type | Location
48  input | - | -
49  output file | - | -
50  output field | vol<Type>Field | $FOAM_CASE/<time>/<outField>
51  \endtable
52 
53  where \c <Type>=Scalar/SymmTensor.
54 
55  References:
56  \verbatim
57  Estimation expressions for L (tag:P), Eq. 10.37:
58  Pope, S. B. (2000).
59  Turbulent flows.
60  Cambridge, UK: Cambridge Univ. Press
61  DOI:10.1017/CBO9780511840531
62  \endverbatim
63 
64 Usage
65  Minimal example by using \c system/controlDict.functions:
66  \verbatim
67  turbulenceFields1
68  {
69  // Mandatory entries (unmodifiable)
70  type turbulenceFields;
71  libs (fieldFunctionObjects);
72 
73  // Mandatory entries (runtime modifiable)
74 
75  // Either of the below
76  // Option-1
77  fields (R devRhoReff);
78 
79  // Option-2
80  field R;
81 
82  // Optional entries (runtime modifiable)
83  prefix <word>;
84 
85  // Inherited entries
86  ...
87  }
88  \endverbatim
89 
90  where the entries mean:
91  \table
92  Property | Description | Type | Reqd | Deflt
93  type | Type name: turbulenceFields | word | yes | -
94  libs | Library name: fieldFunctionObjects | word | yes | -
95  fields | Names of fields to store (see below) | wordList | yes | -
96  field | Name of a field to store (see below) | word | yes | -
97  prefix | Name of output-field prefix | word | no | turbulenceProperties
98  \endtable
99 
100  where \c fields can include:
101  \verbatim
102  k | turbulent kinetic energy
103  epsilon | turbulent kinetic energy dissipation rate
104  omega | specific dissipation rate
105  nuTilda | modified turbulent viscosity
106  nut | turbulent viscosity (incompressible)
107  nuEff | effective turbulent viscosity (incompressible)
108  mut | turbulent viscosity (compressible)
109  muEff | effective turbulent viscosity (compressible)
110  alphat | turbulence thermal diffusivity (compressible)
111  alphaEff | effective turbulence thermal diffusivity (compressible)
112  R | Reynolds stress tensor
113  devReff | deviatoric part of the effective Reynolds stress
114  devRhoReff | divergence of the Reynolds stress
115  L | integral-length scale / mixing-length scale
116  I | turbulence intensity
117  \endverbatim
118 
119  The inherited entries are elaborated in:
120  - \link functionObject.H \endlink
121 
122  Minimal example by using the \c postProcess utility:
123  \verbatim
124  <solver> -postProcess -func turbulenceFields
125  \endverbatim
126 
127 Note
128  - Multiphase applications are not supported.
129  - The governing expression of \c nuTilda is
130  an approximation based on a dimensional analysis.
131 
132 See also
133  - Foam::functionObject
134  - Foam::functionObjects::fvMeshFunctionObject
135  - ExtendedCodeGuide::functionObjects::field::turbulenceFields
136 
137 SourceFiles
138  turbulenceFields.C
139  turbulenceFieldsTemplates.C
140 
141 \*---------------------------------------------------------------------------*/
142 
143 #ifndef functionObjects_turbulenceFields_H
144 #define functionObjects_turbulenceFields_H
145 
146 #include "fvMeshFunctionObject.H"
147 #include "HashSet.H"
148 #include "Enum.H"
149 #include "volFieldsFwd.H"
150 #include "Switch.H"
151 
152 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
153 
154 namespace Foam
155 {
156 namespace functionObjects
157 {
158 
159 /*---------------------------------------------------------------------------*\
160  Class turbulenceFields Declaration
161 \*---------------------------------------------------------------------------*/
162 
163 class turbulenceFields
164 :
165  public fvMeshFunctionObject
166 {
167 public:
168 
169  // Public Enumerations
170 
171  //- Options for the turbulence fields (compressible)
172  enum compressibleField
173  {
174  cfK,
175  cfEpsilon,
176  cfOmega,
177  cfNuTilda,
178  cfMut,
179  cfMuEff,
180  cfAlphat,
181  cfAlphaEff,
182  cfR,
183  cfDevRhoReff,
184  cfL,
185  cfI,
186  cfLESRegion,
187  cffd
188  };
189 
190  //- Names for compressibleField turbulence fields
191  static const Enum<compressibleField> compressibleFieldNames_;
192 
193  //- Options for the turbulence fields (incompressible)
195  {
196  ifK,
197  ifEpsilon,
198  ifOmega,
199  ifNuTilda,
200  ifNut,
201  ifNuEff,
202  ifR,
203  ifDevReff,
204  ifL,
205  ifI,
206  ifLESRegion,
207  iffd
208  };
209 
210  //- Names for incompressibleField turbulence fields
212 
213  //- Name of the turbulence properties dictionary
214  static const word modelName_;
215 
216 
217 protected:
218 
219  // Protected Data
220 
221  //- Flag to track initialisation
222  bool initialised_;
224  //- Name of output-field prefix
227  //- Fields to load
231  // Protected Member Functions
233  //- Unset duplicate fields already registered by other function objects
234  void initialise();
236  //- Return true if compressible turbulence model is identified
237  bool compressible();
238 
239  //- Process the turbulence field
240  template<class Type>
241  void processField
242  (
243  const word& fieldName,
245  );
246 
247  //- Return nuTilda calculated from k and omega
248  template<class Model>
249  tmp<volScalarField> nuTilda(const Model& model) const;
251  //- Return integral length scale, L, calculated from k and epsilon
252  template<class Model>
253  tmp<volScalarField> L(const Model& model) const;
255  //- Return turbulence intensity, I, calculated from k and U
256  template<class Model>
257  tmp<volScalarField> I(const Model& model) const;
260 public:
261 
262  //- Runtime type information
263  TypeName("turbulenceFields");
264 
265 
266  // Constructors
267 
268  //- Construct from Time and dictionary
270  (
271  const word& name,
272  const Time& runTime,
273  const dictionary& dict
274  );
275 
276  //- No copy construct
277  turbulenceFields(const turbulenceFields&) = delete;
278 
279  //- No copy assignment
280  void operator=(const turbulenceFields&) = delete;
282 
283  //- Destructor
284  virtual ~turbulenceFields() = default;
285 
287  // Member Functions
288 
289  //- Read the controls
290  virtual bool read(const dictionary&);
292  //- Calculate turbulence fields
293  virtual bool execute();
294 
295  //- Do nothing.
296  // The turbulence fields are registered and written automatically
297  virtual bool write();
298 };
299 
300 
301 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
302 
303 } // End namespace functionObjects
304 } // End namespace Foam
305 
306 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
307 
308 #ifdef NoRepository
309  #include "turbulenceFieldsTemplates.C"
310 #endif
311 
312 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
313 
314 #endif
315 
316 // ************************************************************************* //
tmp< volScalarField > I(const Model &model) const
Return turbulence intensity, I, calculated from k and U.
dictionary dict
"Effective turbulent dynamic viscosity"
void initialise()
Unset duplicate fields already registered by other function objects.
bool initialised_
Flag to track initialisation.
"Deviatoric part of the effective Reynolds stress"
"Turbulent kinetic energy dissipation rate"
compressibleField
Options for the turbulence fields (compressible)
Forwards and collection of common volume field types.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
turbulenceFields(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
"Effective turbulence thermal diffusivity"
engineTime & runTime
static const word modelName_
Name of the turbulence properties dictionary.
"DES model LES region indicator field"
TypeName("turbulenceFields")
Runtime type information.
"DES model LES region indicator field"
Generic GeometricField class.
Definition: areaFieldsFwd.H:50
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
"Turbulent kinetic energy dissipation rate"
const word & name() const noexcept
Return the name of this functionObject.
wordHashSet fieldSet_
Fields to load.
bool compressible()
Return true if compressible turbulence model is identified.
"Integral-length/Mixing-length scale"
A class for handling words, derived from Foam::string.
Definition: word.H:63
void operator=(const turbulenceFields &)=delete
No copy assignment.
virtual bool read(const dictionary &)
Read the controls.
virtual ~turbulenceFields()=default
Destructor.
static const Enum< incompressibleField > incompressibleFieldNames_
Names for incompressibleField turbulence fields.
static const Enum< compressibleField > compressibleFieldNames_
Names for compressibleField turbulence fields.
Computes various turbulence-related quantities that are not typically output during calculations...
tmp< volScalarField > nuTilda(const Model &model) const
Return nuTilda calculated from k and omega.
void processField(const word &fieldName, const tmp< GeometricField< Type, fvPatchField, volMesh >> &tvalue)
Process the turbulence field.
virtual bool execute()
Calculate turbulence fields.
"Integral-length/Mixing-length scale"
incompressibleField
Options for the turbulence fields (incompressible)
tmp< volScalarField > L(const Model &model) const
Return integral length scale, L, calculated from k and epsilon.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Namespace for OpenFOAM.