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-2021 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 {
37  defineTypeNameAndDebug(fieldAverage, 0);
38  addToRunTimeSelectionTable(functionObject, fieldAverage, dictionary);
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  forAll(faItems_, i)
62  {
63  if (faItems_[i].mean())
64  {
65  if (obr_.found(faItems_[i].meanFieldName()))
66  {
67  obr_.checkOut(*obr_[faItems_[i].meanFieldName()]);
68  }
69  }
70 
71  if (faItems_[i].prime2Mean())
72  {
73  if (obr_.found(faItems_[i].prime2MeanFieldName()))
74  {
75  obr_.checkOut(*obr_[faItems_[i].prime2MeanFieldName()]);
76  }
77  }
78  }
79 }
80 
81 
83 {
84  if (!totalIter_.size())
85  {
86  totalIter_.setSize(faItems_.size(), 1);
87  }
88 
89  if (!totalTime_.size())
90  {
91  totalTime_.setSize(faItems_.size(), obr_.time().deltaTValue());
92  }
93  else
94  {
95  // Check if totalTime_ has been set otherwise initialise
96  forAll(totalTime_, fieldi)
97  {
98  if (totalTime_[fieldi] < 0)
99  {
100  totalTime_[fieldi] = obr_.time().deltaTValue();
101  }
102  }
103  }
104 
105  resetFields();
106 
107  // Add mean fields to the field lists
108  forAll(faItems_, fieldi)
109  {
110  addMeanField<scalar>(fieldi);
111  addMeanField<vector>(fieldi);
112  addMeanField<sphericalTensor>(fieldi);
113  addMeanField<symmTensor>(fieldi);
114  addMeanField<tensor>(fieldi);
115  }
116 
117  // Add prime-squared mean fields to the field lists
118  forAll(faItems_, fieldi)
119  {
120  addPrime2MeanField<scalar, scalar>(fieldi);
121  addPrime2MeanField<vector, symmTensor>(fieldi);
122  }
123 
124  // ensure first averaging works unconditionally
125  prevTimeIndex_ = -1;
126 
127  initialised_ = true;
128 }
129 
130 
132 {
133  Log << " Restarting averaging at time " << obr_.time().timeName()
134  << nl << endl;
135 
136  totalIter_.clear();
137  totalTime_.clear();
138 
139  initialise();
140 }
141 
142 
144 {
145  if (!initialised_)
146  {
147  initialise();
148  }
149 
150  const label currentTimeIndex = obr_.time().timeIndex();
151  const scalar currentTime = obr_.time().value();
152 
153  if (prevTimeIndex_ == currentTimeIndex)
154  {
155  return;
156  }
157  else
158  {
159  prevTimeIndex_ = currentTimeIndex;
160  }
161 
162  if (periodicRestart_ && currentTime > restartPeriod_*periodIndex_)
163  {
164  restart();
165  periodIndex_++;
166  }
167 
168  Log << type() << " " << name() << nl
169  << " Calculating averages" << nl;
170 
171  addMeanSqrToPrime2Mean<scalar, scalar>();
172  addMeanSqrToPrime2Mean<vector, symmTensor>();
173 
174  calculateMeanFields<scalar>();
175  calculateMeanFields<vector>();
176  calculateMeanFields<sphericalTensor>();
177  calculateMeanFields<symmTensor>();
178  calculateMeanFields<tensor>();
179 
180  calculatePrime2MeanFields<scalar, scalar>();
181  calculatePrime2MeanFields<vector, symmTensor>();
182 
183  forAll(faItems_, fieldi)
184  {
185  totalIter_[fieldi]++;
186  totalTime_[fieldi] += obr_.time().deltaTValue();
187  }
188 
189  Log << endl;
190 }
191 
192 
194 {
195  Log << " Writing average fields" << endl;
196 
197  writeFields<scalar>();
198  writeFields<vector>();
199  writeFields<sphericalTensor>();
200  writeFields<symmTensor>();
201  writeFields<tensor>();
202 
203  Log << endl;
204 }
205 
206 
208 {
210  (
211  IOobject
212  (
213  name() + "Properties",
214  obr_.time().timeName(),
215  "uniform",
216  obr_,
219  false
220  )
221  );
222 
223  forAll(faItems_, fieldi)
224  {
225  const word& fieldName = faItems_[fieldi].fieldName();
226  propsDict.add(fieldName, dictionary());
227  propsDict.subDict(fieldName).add("totalIter", totalIter_[fieldi]);
228  propsDict.subDict(fieldName).add("totalTime", totalTime_[fieldi]);
229  }
230 
231  propsDict.regIOobject::write();
232 
233  Log << endl;
234 }
235 
236 
238 {
239  if ((restartOnRestart_ || restartOnOutput_) && log)
240  {
241  Info<< " Starting averaging at time " << obr_.time().timeName()
242  << nl;
243  }
244  else
245  {
246  typeIOobject<timeIOdictionary> propsDictHeader
247  (
248  name() + "Properties",
249  obr_.time().timeName(obr_.time().startTime().value()),
250  "uniform",
251  obr_,
254  false
255  );
256 
257  if (!propsDictHeader.headerOk())
258  {
259  Log << " Starting averaging at time "
260  << obr_.time().timeName() << nl;
261 
262  return;
263  }
264 
265  timeIOdictionary propsDict(propsDictHeader);
266 
267  Log << " Restarting averaging for fields:" << nl;
268 
269  totalIter_.setSize(faItems_.size(), 1);
270 
271  // Initialise totalTime with negative values
272  // to indicate that it has not been set
273  totalTime_.setSize(faItems_.size(), -1);
274 
275  forAll(faItems_, fieldi)
276  {
277  const word& fieldName = faItems_[fieldi].fieldName();
278  if (propsDict.found(fieldName))
279  {
280  dictionary fieldDict(propsDict.subDict(fieldName));
281 
282  totalIter_[fieldi] = fieldDict.lookup<label>("totalIter");
283  totalTime_[fieldi] = fieldDict.lookup<scalar>("totalTime");
284 
285  Log << " " << fieldName
286  << " iters = " << totalIter_[fieldi]
287  << " time = " << totalTime_[fieldi] << nl;
288  }
289  }
290  }
291 }
292 
293 
294 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
295 
297 (
298  const word& name,
299  const Time& runTime,
300  const dictionary& dict
301 )
302 :
303  fvMeshFunctionObject(name, runTime, dict),
304  prevTimeIndex_(-1),
305  restartOnRestart_(false),
306  restartOnOutput_(false),
307  periodicRestart_(false),
308  restartPeriod_(great),
309  initialised_(false),
310  base_(baseType::iter),
311  window_(-1.0),
312  windowName_(""),
313  faItems_(),
314  totalIter_(),
315  totalTime_(),
316  periodIndex_(1)
317 {
318  read(dict);
319 }
320 
321 
322 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
323 
325 {}
326 
327 
328 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
329 
331 {
333 
334  initialised_ = false;
335 
336  Log << type() << " " << name() << ":" << nl;
337 
338  dict.readIfPresent("restartOnRestart", restartOnRestart_);
339  dict.readIfPresent("restartOnOutput", restartOnOutput_);
340  dict.readIfPresent("periodicRestart", periodicRestart_);
341 
342  mean_ = dict.lookupOrDefault<Switch>("mean", true);
343  prime2Mean_ = dict.lookupOrDefault<Switch>("prime2Mean", false);
344  base_ = baseTypeNames_
345  [
346  dict.lookupOrDefault<word>("base", "time")
347  ];
348  window_ = dict.lookupOrDefault<scalar>("window", -1);
349  windowName_ = dict.lookupOrDefault<word>("windowName", "");
350 
351  faItems_ = PtrList<fieldAverageItem>
352  (
353  dict.lookup("fields"),
355  );
356 
357  if (periodicRestart_)
358  {
359  dict.lookup("restartPeriod") >> restartPeriod_;
360  }
361 
362  readAveragingProperties();
363 
364  Log << endl;
365 
366  return true;
367 }
368 
369 
371 {
372  wordList fields(faItems_.size());
373 
374  forAll(faItems_, fieldi)
375  {
376  fields[fieldi] = faItems_[fieldi].fieldName();
377  }
378 
379  return fields;
380 }
381 
382 
384 {
385  calcAverages();
386 
387  return true;
388 }
389 
390 
392 {
393  writeAverages();
394  writeAveragingProperties();
395 
396  if (restartOnOutput_)
397  {
398  restart();
399  }
400 
401  return true;
402 }
403 
404 
405 // ************************************************************************* //
void writeAveragingProperties() const
Write averaging properties - steps and time.
Definition: fieldAverage.C:207
Class used for the read-construction of.
IOdictionary propsDict(systemDict("particleTracksDict", args, runTime))
void initialise()
Reset lists (clear existing values) and initialise averaging.
Definition: fieldAverage.C:82
Calculates the natural logarithm of the specified scalar field.
Definition: log.H:103
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:663
virtual bool read(const dictionary &)
Read the field average data.
Definition: fieldAverage.C:330
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
bool headerOk()
Read header (uses typeGlobalFile to find file) and check.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
Templated form of IOobject providing type information for file reading and header type checking...
Definition: IOobject.H:537
virtual bool write()
Write the field averages.
Definition: fieldAverage.C:391
addToRunTimeSelectionTable(functionObject, Qdot, dictionary)
static const NamedEnum< baseType, 2 > baseTypeNames_
Averaging base type names.
Definition: fieldAverage.H:198
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
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 bool execute()
Calculate the field averages.
Definition: fieldAverage.C:383
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:69
Macros for easy insertion into run-time selection tables.
bool add(entry *, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:1153
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:1002
virtual bool read(const dictionary &)
Read optional controls.
bool read(const char *, int32_t &)
Definition: int32IO.C:85
fieldAverage(const word &name, const Time &runTime, const dictionary &)
Construct from Time and dictionary.
Definition: fieldAverage.C:297
void readAveragingProperties()
Read averaging properties - steps and time.
Definition: fieldAverage.C:237
A class for handling words, derived from string.
Definition: word.H:59
virtual void calcAverages()
Main calculation routine.
Definition: fieldAverage.C:143
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
static const char nl
Definition: Ostream.H:260
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
defineTypeNameAndDebug(Qdot, 0)
timeIOdictionary derived from IOdictionary with globalFile set false to enable writing to processor t...
virtual wordList fields() const
Return the list of fields required.
Definition: fieldAverage.C:370
Info<< "Reading field p_rgh\"<< endl;volScalarField p_rgh(IOobject("p_rgh", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);pressureReference pressureReference(p, p_rgh, pimple.dict(), thermo.incompressible());mesh.schemes().setFluxRequired(p_rgh.name());hydrostaticInitialisation(p_rgh, p, rho, U, gh, ghf, pRef, thermo, pimple.dict());Info<< "Creating field dpdt\"<< endl;volScalarField dpdt(IOobject("dpdt", runTime.timeName(), mesh), mesh, dimensionedScalar(p.dimensions()/dimTime, 0));Info<< "Creating field kinetic energy K\"<< endl;volScalarField K("K", 0.5 *magSqr(U));dimensionedScalar initialMass=fvc::domainIntegrate(rho);multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:131
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
void resetFields()
Checkout fields (causes deletion) from the database.
Definition: fieldAverage.C:59
void restart()
Restart averaging for restartOnOutput.
Definition: fieldAverage.C:131
#define Log
Report write to Foam::Info if the local log switch is true.
messageStream Info
Specialisation of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
virtual ~fieldAverage()
Destructor.
Definition: fieldAverage.C:324
virtual void writeAverages() const
Write averages.
Definition: fieldAverage.C:193
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98
baseType
Enumeration defining the averaging base type.
Definition: fieldAverage.H:168
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:864