fieldAverage.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 | Copyright (C) 2011-2016 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 "volFields.H"
28 #include "fieldAverageItem.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 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
44 
46 {
47  forAll(faItems_, i)
48  {
49  if (faItems_[i].mean())
50  {
51  if (obr_.found(faItems_[i].meanFieldName()))
52  {
53  obr_.checkOut(*obr_[faItems_[i].meanFieldName()]);
54  }
55  }
56 
57  if (faItems_[i].prime2Mean())
58  {
59  if (obr_.found(faItems_[i].prime2MeanFieldName()))
60  {
61  obr_.checkOut(*obr_[faItems_[i].prime2MeanFieldName()]);
62  }
63  }
64  }
65 }
66 
67 
69 {
70  resetFields();
71 
72  // Add mean fields to the field lists
73  forAll(faItems_, fieldi)
74  {
75  addMeanField<scalar>(fieldi);
76  addMeanField<vector>(fieldi);
77  addMeanField<sphericalTensor>(fieldi);
78  addMeanField<symmTensor>(fieldi);
79  addMeanField<tensor>(fieldi);
80  }
81 
82  // Add prime-squared mean fields to the field lists
83  forAll(faItems_, fieldi)
84  {
85  addPrime2MeanField<scalar, scalar>(fieldi);
86  addPrime2MeanField<vector, symmTensor>(fieldi);
87  }
88 
89  // ensure first averaging works unconditionally
90  prevTimeIndex_ = -1;
91 
92  initialised_ = true;
93 }
94 
95 
97 {
98  Log
99  << " Restarting averaging at time " << obr_.time().timeName()
100  << nl << endl;
101 
102  totalIter_.clear();
103  totalIter_.setSize(faItems_.size(), 1);
104 
105  totalTime_.clear();
106  totalTime_.setSize(faItems_.size(), obr_.time().deltaTValue());
107 
108  initialize();
109 }
110 
111 
113 {
114  if (!initialised_)
115  {
116  initialize();
117  }
118 
119  const label currentTimeIndex = obr_.time().timeIndex();
120  const scalar currentTime = obr_.time().value();
121 
122  if (prevTimeIndex_ == currentTimeIndex)
123  {
124  return;
125  }
126  else
127  {
128  prevTimeIndex_ = currentTimeIndex;
129  }
130 
131  if (periodicRestart_ && currentTime > restartPeriod_*periodIndex_)
132  {
133  restart();
134  periodIndex_++;
135  }
136 
137  Log
138  << type() << " " << name() << " write:" << nl
139  << " Calculating averages" << nl;
140 
141  addMeanSqrToPrime2Mean<scalar, scalar>();
142  addMeanSqrToPrime2Mean<vector, symmTensor>();
143 
144  calculateMeanFields<scalar>();
145  calculateMeanFields<vector>();
146  calculateMeanFields<sphericalTensor>();
147  calculateMeanFields<symmTensor>();
148  calculateMeanFields<tensor>();
149 
150  calculatePrime2MeanFields<scalar, scalar>();
151  calculatePrime2MeanFields<vector, symmTensor>();
152 
153  forAll(faItems_, fieldi)
154  {
155  totalIter_[fieldi]++;
156  totalTime_[fieldi] += obr_.time().deltaTValue();
157  }
158 
159  Log << endl;
160 }
161 
162 
164 {
165  Log << " Writing average fields" << endl;
166 
167  writeFields<scalar>();
168  writeFields<vector>();
169  writeFields<sphericalTensor>();
170  writeFields<symmTensor>();
171  writeFields<tensor>();
172 
173  Log << endl;
174 }
175 
176 
178 {
180  (
181  IOobject
182  (
183  name() + "Properties",
184  obr_.time().timeName(),
185  "uniform",
186  obr_,
189  false
190  )
191  );
192 
193  forAll(faItems_, fieldi)
194  {
195  const word& fieldName = faItems_[fieldi].fieldName();
196  propsDict.add(fieldName, dictionary());
197  propsDict.subDict(fieldName).add("totalIter", totalIter_[fieldi]);
198  propsDict.subDict(fieldName).add("totalTime", totalTime_[fieldi]);
199  }
200 
201  propsDict.regIOobject::write();
202 
203  Log << endl;
204 }
205 
206 
208 {
209  totalIter_.clear();
210  totalIter_.setSize(faItems_.size(), 1);
211 
212  totalTime_.clear();
213  totalTime_.setSize(faItems_.size(), obr_.time().deltaTValue());
214 
215  if ((restartOnRestart_ || restartOnOutput_) && log)
216  {
217  Info<< " Starting averaging at time " << obr_.time().timeName()
218  << nl;
219  }
220  else
221  {
222  IOobject propsDictHeader
223  (
224  name() + "Properties",
225  obr_.time().timeName(obr_.time().startTime().value()),
226  "uniform",
227  obr_,
230  false
231  );
232 
233  if (!propsDictHeader.headerOk())
234  {
235  Log
236  << " Starting averaging at time "
237  << obr_.time().timeName() << nl;
238 
239  return;
240  }
241 
242  IOdictionary propsDict(propsDictHeader);
243 
244  Log << " Restarting averaging for fields:" << nl;
245 
246  forAll(faItems_, fieldi)
247  {
248  const word& fieldName = faItems_[fieldi].fieldName();
249  if (propsDict.found(fieldName))
250  {
251  dictionary fieldDict(propsDict.subDict(fieldName));
252 
253  totalIter_[fieldi] = readLabel(fieldDict.lookup("totalIter"));
254  totalTime_[fieldi] = readScalar(fieldDict.lookup("totalTime"));
255 
256  Log
257  << " " << fieldName
258  << " iters = " << totalIter_[fieldi]
259  << " time = " << totalTime_[fieldi] << nl;
260  }
261  }
262  }
263 }
264 
265 
266 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
267 
269 (
270  const word& name,
271  const Time& runTime,
272  const dictionary& dict
273 )
274 :
275  fvMeshFunctionObject(name, runTime, dict),
276  prevTimeIndex_(-1),
277  restartOnRestart_(false),
278  restartOnOutput_(false),
279  periodicRestart_(false),
280  restartPeriod_(GREAT),
281  initialised_(false),
282  faItems_(),
283  totalIter_(),
284  totalTime_(),
285  periodIndex_(1)
286 {
287  read(dict);
288 }
289 
290 
291 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
292 
294 {}
295 
296 
297 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
298 
300 {
302 
303  initialised_ = false;
304 
305  Log << type() << " " << name() << ":" << nl;
306 
307  dict.readIfPresent("restartOnRestart", restartOnRestart_);
308  dict.readIfPresent("restartOnOutput", restartOnOutput_);
309  dict.readIfPresent("periodicRestart", periodicRestart_);
310  dict.lookup("fields") >> faItems_;
311 
312  if (periodicRestart_)
313  {
314  dict.lookup("restartPeriod") >> restartPeriod_;
315  }
316 
317  readAveragingProperties();
318 
319  Log << endl;
320 
321  return true;
322 }
323 
324 
326 {
327  calcAverages();
328 
329  return true;
330 }
331 
332 
334 {
335  writeAverages();
336  writeAveragingProperties();
337 
338  if (restartOnOutput_)
339  {
340  restart();
341  }
342 
343  return true;
344 }
345 
346 
347 // ************************************************************************* //
virtual bool read(const dictionary &)
Read the field average data.
Definition: fieldAverage.C:299
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
fieldAverage(const fieldAverage &)
Disallow default bitwise copy construct.
dimensionedScalar log(const dimensionedScalar &ds)
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
virtual bool write()
Write the field averages.
Definition: fieldAverage.C:333
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:633
IOdictionary propsDict(IOobject("particleTrackProperties", runTime.constant(), mesh, IOobject::MUST_READ_IF_MODIFIED))
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
virtual bool execute()
Calculate the field averages.
Definition: fieldAverage.C:325
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
Macros for easy insertion into run-time selection tables.
bool add(entry *, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:737
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
virtual void writeAverages() const
Write averages.
Definition: fieldAverage.C:163
virtual bool read(const dictionary &)
Read optional controls.
bool read(const char *, int32_t &)
Definition: int32IO.C:85
void writeAveragingProperties() const
Write averaging properties - steps and time.
Definition: fieldAverage.C:177
void readAveragingProperties()
Read averaging properties - steps and time.
Definition: fieldAverage.C:207
A class for handling words, derived from string.
Definition: word.H:59
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
virtual void calcAverages()
Main calculation routine.
Definition: fieldAverage.C:112
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:306
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
label readLabel(Istream &is)
Definition: label.H:64
static const char nl
Definition: Ostream.H:262
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
void initialize()
Reset lists (clear existing values) and initialize averaging.
Definition: fieldAverage.C:68
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:461
addToRunTimeSelectionTable(functionObject, blendingFactor, dictionary)
void resetFields()
Checkout fields (causes deletion) from the database.
Definition: fieldAverage.C:45
void restart()
Restart averaging for restartOnOutput.
Definition: fieldAverage.C:96
#define Log
Report write to Foam::Info if the local log switch is true.
messageStream Info
bool headerOk()
Read and check header info.
Definition: IOobject.C:400
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
defineTypeNameAndDebug(fvMeshFunctionObject, 0)
virtual ~fieldAverage()
Destructor.
Definition: fieldAverage.C:293
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:451