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-2015 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 "Time.H"
29 
30 #include "fieldAverageItem.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  defineTypeNameAndDebug(fieldAverage, 0);
37 }
38 
39 
40 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
41 
43 {
44  forAll(faItems_, i)
45  {
46  if (faItems_[i].mean())
47  {
48  if (obr_.found(faItems_[i].meanFieldName()))
49  {
50  obr_.checkOut(*obr_[faItems_[i].meanFieldName()]);
51  }
52  }
53 
54  if (faItems_[i].prime2Mean())
55  {
56  if (obr_.found(faItems_[i].prime2MeanFieldName()))
57  {
58  obr_.checkOut(*obr_[faItems_[i].prime2MeanFieldName()]);
59  }
60  }
61  }
62 }
63 
64 
66 {
67  resetFields();
68 
69  Info<< type() << " " << name_ << ":" << nl;
70 
71  // Add mean fields to the field lists
72  forAll(faItems_, fieldI)
73  {
74  addMeanField<scalar>(fieldI);
75  addMeanField<vector>(fieldI);
76  addMeanField<sphericalTensor>(fieldI);
77  addMeanField<symmTensor>(fieldI);
78  addMeanField<tensor>(fieldI);
79  }
80 
81  // Add prime-squared mean fields to the field lists
82  forAll(faItems_, fieldI)
83  {
84  addPrime2MeanField<scalar, scalar>(fieldI);
85  addPrime2MeanField<vector, symmTensor>(fieldI);
86  }
87 
88  forAll(faItems_, fieldI)
89  {
90  if (!faItems_[fieldI].active())
91  {
92  WarningIn("void Foam::fieldAverage::initialize()")
93  << "Field " << faItems_[fieldI].fieldName()
94  << " not found in database for averaging";
95  }
96  }
97 
98  // ensure first averaging works unconditionally
99  prevTimeIndex_ = -1;
100 
101  Info<< endl;
102 
103  initialised_ = true;
104 }
105 
106 
108 {
109  if (!initialised_)
110  {
111  initialize();
112  }
113 
114  const label currentTimeIndex =
115  static_cast<const fvMesh&>(obr_).time().timeIndex();
116 
117  if (prevTimeIndex_ == currentTimeIndex)
118  {
119  return;
120  }
121  else
122  {
123  prevTimeIndex_ = currentTimeIndex;
124  }
125 
126  Info<< type() << " " << name_ << " output:" << nl;
127 
128  Info<< " Calculating averages" << nl;
129 
130  addMeanSqrToPrime2Mean<scalar, scalar>();
131  addMeanSqrToPrime2Mean<vector, symmTensor>();
132 
133  calculateMeanFields<scalar>();
134  calculateMeanFields<vector>();
135  calculateMeanFields<sphericalTensor>();
136  calculateMeanFields<symmTensor>();
137  calculateMeanFields<tensor>();
138 
139  calculatePrime2MeanFields<scalar, scalar>();
140  calculatePrime2MeanFields<vector, symmTensor>();
141 
142  forAll(faItems_, fieldI)
143  {
144  totalIter_[fieldI]++;
145  totalTime_[fieldI] += obr_.time().deltaTValue();
146  }
147 }
148 
149 
151 {
152  Info<< " Writing average fields" << endl;
153 
154  writeFields<scalar>();
155  writeFields<vector>();
156  writeFields<sphericalTensor>();
157  writeFields<symmTensor>();
158  writeFields<tensor>();
159 }
160 
161 
163 {
165  (
166  IOobject
167  (
168  "fieldAveragingProperties",
169  obr_.time().timeName(),
170  "uniform",
171  obr_,
174  false
175  )
176  );
177 
178  forAll(faItems_, fieldI)
179  {
180  const word& fieldName = faItems_[fieldI].fieldName();
181  propsDict.add(fieldName, dictionary());
182  propsDict.subDict(fieldName).add("totalIter", totalIter_[fieldI]);
183  propsDict.subDict(fieldName).add("totalTime", totalTime_[fieldI]);
184  }
185 
186  propsDict.regIOobject::write();
187 }
188 
189 
191 {
192  totalIter_.clear();
193  totalIter_.setSize(faItems_.size(), 1);
194 
195  totalTime_.clear();
196  totalTime_.setSize(faItems_.size(), obr_.time().deltaTValue());
197 
198  if (resetOnRestart_ || resetOnOutput_)
199  {
200  Info<< " Starting averaging at time " << obr_.time().timeName()
201  << nl;
202  }
203  else
204  {
205  IOobject propsDictHeader
206  (
207  "fieldAveragingProperties",
208  obr_.time().timeName(obr_.time().startTime().value()),
209  "uniform",
210  obr_,
213  false
214  );
215 
216  if (!propsDictHeader.headerOk())
217  {
218  Info<< " Starting averaging at time " << obr_.time().timeName()
219  << nl;
220  return;
221  }
222 
223  IOdictionary propsDict(propsDictHeader);
224 
225  Info<< " Restarting averaging for fields:" << nl;
226  forAll(faItems_, fieldI)
227  {
228  const word& fieldName = faItems_[fieldI].fieldName();
229  if (propsDict.found(fieldName))
230  {
231  dictionary fieldDict(propsDict.subDict(fieldName));
232 
233  totalIter_[fieldI] = readLabel(fieldDict.lookup("totalIter"));
234  totalTime_[fieldI] = readScalar(fieldDict.lookup("totalTime"));
235  Info<< " " << fieldName
236  << " iters = " << totalIter_[fieldI]
237  << " time = " << totalTime_[fieldI] << nl;
238  }
239  }
240  }
241 }
242 
243 
244 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
245 
247 (
248  const word& name,
249  const objectRegistry& obr,
250  const dictionary& dict,
251  const bool loadFromFiles
252 )
253 :
254  name_(name),
255  obr_(obr),
256  active_(true),
257  prevTimeIndex_(-1),
258  resetOnRestart_(false),
259  resetOnOutput_(false),
260  initialised_(false),
261  faItems_(),
262  totalIter_(),
263  totalTime_()
264 {
265  // Only active if a fvMesh is available
266  if (isA<fvMesh>(obr_))
267  {
268  read(dict);
269  }
270  else
271  {
272  active_ = false;
273  WarningIn
274  (
275  "fieldAverage::fieldAverage"
276  "("
277  "const word&, "
278  "const objectRegistry&, "
279  "const dictionary&, "
280  "const bool "
281  ")"
282  ) << "No fvMesh available, deactivating " << name_ << nl
283  << endl;
284  }
285 }
286 
287 
288 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
289 
291 {}
292 
293 
294 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
295 
297 {
298  if (active_)
299  {
300  initialised_ = false;
301 
302  Info<< type() << " " << name_ << ":" << nl;
303 
304  dict.readIfPresent("resetOnRestart", resetOnRestart_);
305  dict.readIfPresent("resetOnOutput", resetOnOutput_);
306  dict.lookup("fields") >> faItems_;
307 
308  readAveragingProperties();
309 
310  Info<< endl;
311  }
312 }
313 
314 
316 {
317  if (active_)
318  {
319  calcAverages();
320  Info<< endl;
321  }
322 }
323 
324 
326 {
327  if (active_)
328  {
329  calcAverages();
330  Info<< endl;
331  }
332 }
333 
334 
336 {}
337 
338 
340 {
341  if (active_)
342  {
343  writeAverages();
344  writeAveragingProperties();
345 
346  if (resetOnOutput_)
347  {
348  Info<< " Restarting averaging at time " << obr_.time().timeName()
349  << nl << endl;
350 
351  totalIter_.clear();
352  totalIter_.setSize(faItems_.size(), 1);
353 
354  totalTime_.clear();
355  totalTime_.setSize(faItems_.size(), obr_.time().deltaTValue());
356 
357  initialize();
358  }
359 
360  Info<< endl;
361  }
362 }
363 
364 
366 {
367  // Do nothing
368 }
369 
370 
372 {
373  // Do nothing
374 }
375 
376 
377 // ************************************************************************* //
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:306
fieldAverage(const fieldAverage &)
Disallow default bitwise copy construct.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
virtual ~fieldAverage()
Destructor.
Definition: fieldAverage.C:290
virtual void calcAverages()
Main calculation routine.
Definition: fieldAverage.C:107
A class for handling words, derived from string.
Definition: word.H:59
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
messageStream Info
virtual void timeSet()
Called when time was set at the end of the Time::operator++.
Definition: fieldAverage.C:335
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:638
Namespace for OpenFOAM.
label readLabel(Istream &is)
Definition: label.H:64
void writeAveragingProperties() const
Write averaging properties - steps and time.
Definition: fieldAverage.C:162
static const char nl
Definition: Ostream.H:260
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
#define WarningIn(functionName)
Report a warning using Foam::Warning.
IOdictionary propsDict(IOobject( "particleTrackProperties", runTime.constant(), mesh, IOobject::MUST_READ_IF_MODIFIED ))
void readAveragingProperties()
Read averaging properties - steps and time.
Definition: fieldAverage.C:190
#define forAll(list, i)
Definition: UList.H:421
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
void resetFields()
Checkout fields (causes deletion) from the database.
Definition: fieldAverage.C:42
virtual void write()
Calculate the field average data and write.
Definition: fieldAverage.C:339
bool add(entry *, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:739
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:589
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:452
Registry of regIOobjects.
void initialize()
Reset lists (clear existing values) and initialize averaging.
Definition: fieldAverage.C:65
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
bool read(const char *, int32_t &)
Definition: int32IO.C:87
virtual void updateMesh(const mapPolyMesh &)
Update mesh.
Definition: fieldAverage.C:365
bool headerOk()
Read and check header info.
Definition: IOobject.C:424
virtual void end()
Execute the averaging at the final time-loop, currently does nothing.
Definition: fieldAverage.C:325
virtual void writeAverages() const
Write averages.
Definition: fieldAverage.C:150
defineTypeNameAndDebug(combustionModel, 0)
virtual void movePoints(const polyMesh &)
Move points.
Definition: fieldAverage.C:371
virtual void read(const dictionary &)
Read the field average data.
Definition: fieldAverage.C:296
virtual void execute()
Execute the averaging.
Definition: fieldAverage.C:315