fieldAverage.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) 2011-2023 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 "fieldAverage.H"
27 #include "fieldAverageItem.H"
28 #include "timeIOdictionary.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace functionObjects
36 {
39 }
40 }
41 
42 
43 template<>
44 const char* Foam::NamedEnum
45 <
47  2
48 >::names[] = { "iteration", "time"};
49 
50 const Foam::NamedEnum
51 <
53  2
55 
56 
57 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
58 
60 {
61  // Initialise any unset times
62  forAll(totalTime_, fieldi)
63  {
64  if (totalTime_[fieldi] < 0)
65  {
66  totalTime_[fieldi] = 0;
67  }
68  }
69 
70  // Initialise mean fields
71  forAll(faItems_, fieldi)
72  {
73  initialiseMeanField<scalar>(fieldi);
74  initialiseMeanField<vector>(fieldi);
75  initialiseMeanField<sphericalTensor>(fieldi);
76  initialiseMeanField<symmTensor>(fieldi);
77  initialiseMeanField<tensor>(fieldi);
78  }
79  forAll(faItems_, fieldi)
80  {
81  initialisePrime2MeanField<scalar, scalar>(fieldi);
82  initialisePrime2MeanField<vector, symmTensor>(fieldi);
83  }
84 }
85 
86 
88 {
89  Log << " Restarting averaging at time " << time_.name() << nl;
90 
91  // Clear the times
92  totalIter_ = 0;
93  totalTime_ = -1;
94 
95  // Clear mean fields
96  forAll(faItems_, i)
97  {
98  if (faItems_[i].mean())
99  {
100  if (obr_.found(faItems_[i].meanFieldName()))
101  {
102  obr_.checkOut(*obr_[faItems_[i].meanFieldName()]);
103  }
104  }
105 
106  if (faItems_[i].prime2Mean())
107  {
108  if (obr_.found(faItems_[i].prime2MeanFieldName()))
109  {
110  obr_.checkOut(*obr_[faItems_[i].prime2MeanFieldName()]);
111  }
112  }
113  }
114 
115  // Re-create any mean fields
116  initialise();
117 
118  // Ensure first averaging works unconditionally
119  prevTimeIndex_ = -1;
120 }
121 
122 
124 {
125  Log << type() << " " << name() << ":" << nl;
126 
127  const label currentTimeIndex = time_.timeIndex();
128  const scalar currentTime = time_.value();
129 
130  if (prevTimeIndex_ == currentTimeIndex)
131  {
132  return;
133  }
134 
135  prevTimeIndex_ = currentTimeIndex;
136 
137  if (periodicRestart_ && currentTime > restartPeriod_*periodIndex_)
138  {
139  restart();
140  periodIndex_++;
141  }
142  else
143  {
144  initialise();
145  }
146 
147  // Increment the time and iteration totals
148  forAll(faItems_, fieldi)
149  {
150  totalIter_[fieldi]++;
151  totalTime_[fieldi] += time_.deltaTValue();
152  }
153 
154  Log << " Calculating averages" << nl;
155 
156  addMeanSqrToPrime2Mean<scalar, scalar>();
157  addMeanSqrToPrime2Mean<vector, symmTensor>();
158 
159  calculateMeanFields<scalar>();
160  calculateMeanFields<vector>();
161  calculateMeanFields<sphericalTensor>();
162  calculateMeanFields<symmTensor>();
163  calculateMeanFields<tensor>();
164 
165  calculatePrime2MeanFields<scalar, scalar>();
166  calculatePrime2MeanFields<vector, symmTensor>();
167 
168  Log << endl;
169 }
170 
171 
173 {
174  Log << type() << " " << name() << ":" << nl
175  << " Writing average fields" << endl;
176 
177  writeFields<scalar>();
178  writeFields<vector>();
179  writeFields<sphericalTensor>();
180  writeFields<symmTensor>();
181  writeFields<tensor>();
182 
184  (
185  IOobject
186  (
187  name() + "Properties",
188  time_.name(),
189  "uniform",
190  obr_,
193  false
194  )
195  );
196 
197  forAll(faItems_, fieldi)
198  {
199  const word& fieldName = faItems_[fieldi].fieldName();
200  propsDict.add(fieldName, dictionary());
201  propsDict.subDict(fieldName).add("totalIter", totalIter_[fieldi]);
202  propsDict.subDict(fieldName).add("totalTime", totalTime_[fieldi]);
203  }
204 
205  propsDict.regIOobject::write();
206 
207  Log << endl;
208 }
209 
210 
212 (
213  const dictionary& dict,
214  const bool construct
215 )
216 {
217  dict.readIfPresent("restartOnRestart", restartOnRestart_);
218  dict.readIfPresent("restartOnOutput", restartOnOutput_);
219  dict.readIfPresent("periodicRestart", periodicRestart_);
220 
221  if (periodicRestart_)
222  {
223  dict.lookup("restartPeriod") >> restartPeriod_;
224  }
225 
226  mean_ = dict.lookupOrDefault<Switch>("mean", true);
227  prime2Mean_ = dict.lookupOrDefault<Switch>("prime2Mean", false);
228  base_ = baseTypeNames_[dict.lookupOrDefault<word>("base", "time")];
229  window_ = dict.lookupOrDefault<scalar>("window", -1);
230  windowName_ = dict.lookupOrDefault<word>("windowName", "");
231 
232  if (construct)
233  {
234  // First read of a run. Look for properties dict and read total
235  // iter/time values from it if available.
236 
237  faItems_.clear();
238  totalIter_.clear();
239  totalTime_.clear();
240 
241  faItems_ =
243  (
244  dict.lookup("fields"),
246  );
247  totalIter_.setSize(faItems_.size(), 0);
248  totalTime_.setSize(faItems_.size(), -1);
249 
250  typeIOobject<timeIOdictionary> propsDictHeader
251  (
252  name() + "Properties",
253  time_.name(),
254  "uniform",
255  obr_,
258  false
259  );
260 
262  if
263  (
264  !restartOnRestart_
265  && !restartOnOutput_
266  && propsDictHeader.headerOk()
267  )
268  {
269  propsDict = timeIOdictionary(propsDictHeader);
270  }
271 
272  bool first = true;
273  forAll(faItems_, fieldi)
274  {
275  const word& fieldName = faItems_[fieldi].fieldName();
276  if (!propsDict.found(fieldName))
277  {
278  if (first)
279  {
280  Log << " Starting averaging for fields:" << nl;
281  first = false;
282  }
283  Log << " " << fieldName << nl;
284  }
285  }
286 
287  first = true;
288  forAll(faItems_, fieldi)
289  {
290  const word& fieldName = faItems_[fieldi].fieldName();
291  if (propsDict.found(fieldName))
292  {
293  if (first)
294  {
295  Log << " Restarting averaging for fields:" << nl;
296  first = false;
297  }
298  const dictionary& fieldPropsDict = propsDict.subDict(fieldName);
299  totalIter_[fieldi] = fieldPropsDict.lookup<label>("totalIter");
300  totalTime_[fieldi] = fieldPropsDict.lookup<scalar>("totalTime");
301  Log << " " << fieldName
302  << " iters = " << totalIter_[fieldi]
303  << " time = " << totalTime_[fieldi] << nl;
304  }
305  }
306  }
307  else
308  {
309  // Re-read during a run. Read the items and copy the per-field total
310  // iter/times from the old items to the new.
311 
312  PtrList<fieldAverageItem> faItems0;
313  List<label> totalIter0;
314  List<scalar> totalTime0;
315  faItems0.transfer(faItems_);
316  totalIter0.transfer(totalIter_);
317  totalTime0.transfer(totalTime_);
318 
319  faItems_ =
321  (
322  dict.lookup("fields"),
324  );
325  totalIter_.resize(faItems_.size(), 0);
326  totalTime_.resize(faItems_.size(), -1);
327 
328  // Map from field to old-field index
329  labelList fieldiFieldi0s(faItems_.size(), -1);
330  forAll(faItems_, fieldi)
331  {
332  const word& fieldName = faItems_[fieldi].fieldName();
333  forAll(faItems0, fieldi0)
334  {
335  if (faItems0[fieldi0].fieldName() == fieldName)
336  {
337  fieldiFieldi0s[fieldi] = fieldi0;
338  break;
339  }
340  }
341  }
342 
343  bool first = true;
344  forAll(faItems_, fieldi)
345  {
346  if (fieldiFieldi0s[fieldi] == -1)
347  {
348  if (first)
349  {
350  Log << " Starting averaging for fields:" << nl;
351  first = true;
352  }
353  Log << " " << faItems_[fieldi].fieldName() << nl;
354  }
355  }
356 
357  first = true;
358  forAll(faItems_, fieldi)
359  {
360  if (fieldiFieldi0s[fieldi] != -1)
361  {
362  if (first)
363  {
364  Log << " Continuing averaging for fields:" << nl;
365  first = true;
366  }
367  totalIter_[fieldi] = totalIter0[fieldiFieldi0s[fieldi]];
368  totalTime_[fieldi] = totalTime0[fieldiFieldi0s[fieldi]];
369  Log << " " << faItems_[fieldi].fieldName()
370  << " iters = " << totalIter_[fieldi]
371  << " time = " << totalTime_[fieldi] << nl;
372  }
373  }
374  }
375 }
376 
377 
378 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
379 
381 (
382  const word& name,
383  const Time& runTime,
384  const dictionary& dict
385 )
386 :
387  fvMeshFunctionObject(name, runTime, dict),
388  prevTimeIndex_(-1),
389  restartOnRestart_(false),
390  restartOnOutput_(false),
391  periodicRestart_(false),
392  restartPeriod_(great),
393  base_(baseType::iter),
394  window_(-1.0),
395  windowName_(""),
396  faItems_(),
397  totalIter_(),
398  totalTime_(),
399  periodIndex_(1)
400 {
402 
403  Log << type() << " " << name << ":" << nl;
404 
405  read(dict, true);
406 
407  // Read any available mean fields
408  forAll(faItems_, fieldi)
409  {
410  readMeanField<scalar>(fieldi);
411  readMeanField<vector>(fieldi);
412  readMeanField<sphericalTensor>(fieldi);
413  readMeanField<symmTensor>(fieldi);
414  readMeanField<tensor>(fieldi);
415  }
416  forAll(faItems_, fieldi)
417  {
418  readPrime2MeanField<scalar, scalar>(fieldi);
419  readPrime2MeanField<vector, symmTensor>(fieldi);
420  }
421 
422  Log << endl;
423 }
424 
425 
426 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
427 
429 {}
430 
431 
432 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
433 
435 {
437 
438  Log << type() << " " << name() << ":" << nl;
439 
440  read(dict, false);
441 
442  Log << endl;
443 
444  return true;
445 }
446 
447 
449 {
450  wordList fields(faItems_.size());
451 
452  forAll(faItems_, fieldi)
453  {
454  fields[fieldi] = faItems_[fieldi].fieldName();
455  }
456 
457  return fields;
458 }
459 
460 
462 {
464  {
466  << "fieldAverage is not supported with the foamPostProcess utility"
467  << endl;
468 
469  return false;
470  }
471 
472  calcAverages();
473 
474  return true;
475 }
476 
477 
479 {
480  writeAverages();
481 
482  if (restartOnOutput_)
483  {
484  restart();
485  }
486 
487  return true;
488 }
489 
490 
491 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
@ MUST_READ_IF_MODIFIED
Definition: IOobject.H:119
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Definition: List.C:342
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:54
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
void transfer(PtrList< T > &)
Transfer the contents of the argument PtrList into this PtrList.
Definition: PtrList.C:189
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
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:860
Abstract base-class for Time/database functionObjects.
virtual const word & type() const =0
Runtime type information.
static bool postProcess
Global post-processing mode switch.
const word & name() const
Return the name of this functionObject.
Class used for the read-construction of.
Calculates average quantities for a user-specified selection of volumetric and surface fields.
Definition: fieldAverage.H:163
baseType
Enumeration defining the averaging base type.
Definition: fieldAverage.H:168
void initialise()
Reset lists (clear existing values) and initialise averaging.
Definition: fieldAverage.C:59
void restart()
Restart averaging for restartOnOutput.
Definition: fieldAverage.C:87
static const NamedEnum< baseType, 2 > baseTypeNames_
Averaging base type names.
Definition: fieldAverage.H:194
virtual void calcAverages()
Main calculation routine.
Definition: fieldAverage.C:123
virtual wordList fields() const
Return the list of fields required.
Definition: fieldAverage.C:448
virtual ~fieldAverage()
Destructor.
Definition: fieldAverage.C:428
void read(const dictionary &dict, const bool construct)
Read.
Definition: fieldAverage.C:212
List< scalar > totalTime_
Total time counter.
Definition: fieldAverage.H:219
PtrList< fieldAverageItem > faItems_
List of field average items, describing what averages to be.
Definition: fieldAverage.H:213
virtual void writeAverages() const
Write averages.
Definition: fieldAverage.C:172
virtual bool execute()
Calculate the field averages.
Definition: fieldAverage.C:461
virtual bool write()
Write the field averages.
Definition: fieldAverage.C:478
fieldAverage(const word &name, const Time &runTime, const dictionary &)
Construct from Time and dictionary.
Definition: fieldAverage.C:381
Specialisation of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
virtual bool read(const dictionary &)
Read optional controls.
timeIOdictionary derived from IOdictionary with globalFile set false to enable writing to processor t...
Templated form of IOobject providing type information for file reading and header type checking.
Definition: IOobject.H:531
bool headerOk()
Read header (uses typeGlobalFile to find file) and check.
A class for handling words, derived from string.
Definition: word.H:62
Info<< "Calculating turbulent flame speed field St\n"<< endl;volScalarField St(IOobject("St", runTime.name(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), flameWrinkling->Xi() *Su);multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:230
#define WarningInFunction
Report a warning using Foam::Warning.
#define Log
Report write to Foam::Info if the local log switch is true.
defineTypeNameAndDebug(adjustTimeStepToCombustion, 0)
addToRunTimeSelectionTable(functionObject, adjustTimeStepToCombustion, dictionary)
Namespace for OpenFOAM.
bool read(const char *, int32_t &)
Definition: int32IO.C:85
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
labelList first(const UList< labelPair > &p)
Definition: patchToPatch.C:39
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
static const char nl
Definition: Ostream.H:260
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
dictionary dict
IOdictionary propsDict(systemDict("particleTracksDict", args, runTime))