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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2015-2017 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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\*---------------------------------------------------------------------------*/
28
29#include "fieldAverage.H"
30#include "volFields.H"
31#include "fieldAverageItem.H"
33
34// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35
36namespace Foam
37{
38namespace functionObjects
39{
42}
43}
44
45
46// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
47
49{
50 for (fieldAverageItem& item : faItems_)
51 {
52 // Note: not clearing data needed for restart
53 item.clear(obr(), false);
54 }
55
56 Log << type() << " " << name() << ":" << nl;
57
58 // Add mean fields to the field lists
59 for (fieldAverageItem& item : faItems_)
60 {
61 addMeanField<scalar>(item);
62 addMeanField<vector>(item);
63 addMeanField<sphericalTensor>(item);
64 addMeanField<symmTensor>(item);
65 addMeanField<tensor>(item);
66 }
67
68 // Add prime-squared mean fields to the field lists
69 for (fieldAverageItem& item : faItems_)
70 {
71 addPrime2MeanField<scalar, scalar>(item);
72 addPrime2MeanField<vector, symmTensor>(item);
73 }
74
75 // Add window fields to the field lists
76 for (const fieldAverageItem& item : faItems_)
77 {
78 restoreWindowFields<scalar>(item);
79 restoreWindowFields<vector>(item);
80 restoreWindowFields<sphericalTensor>(item);
81 restoreWindowFields<symmTensor>(item);
82 restoreWindowFields<tensor>(item);
83 }
84
85
86 for (const fieldAverageItem& item : faItems_)
87 {
88 if (!item.active())
89 {
91 << "Field " << item.fieldName()
92 << " not found in database for averaging";
93 }
94 }
95
96 // Ensure first averaging works unconditionally
97 prevTimeIndex_ = -1;
98
99 Log << endl;
100 initialised_ = true;
101}
102
103
105{
106 Log << " Restarting averaging at time "
107 << obr().time().timeOutputValue()
108 << nl << endl;
109
110 for (fieldAverageItem& item : faItems_)
111 {
112 item.clear(obr(), true);
113 }
114
115 initialize();
116}
117
118
120{
121 if (!initialised_)
122 {
123 initialize();
124 }
125
126 const label currentTimeIndex = obr().time().timeIndex();
127 const scalar currentTime = obr().time().value();
128
129 if (prevTimeIndex_ == currentTimeIndex)
130 {
131 return;
132 }
133 else
134 {
135 prevTimeIndex_ = currentTimeIndex;
136 }
137
138 bool doRestart = false;
139 if (periodicRestart_)
140 {
141 const scalar deltaT = obr().time().deltaTValue();
142 const scalar nextTime = restartPeriod_*periodIndex_ + 0.5*deltaT;
143
144 if (currentTime > nextTime)
145 {
146 doRestart = true;
147 ++periodIndex_;
148 }
149 }
150
151 if (currentTime >= restartTime_)
152 {
153 doRestart = true; // Restart is overdue.
154 restartTime_ = GREAT; // Avoid triggering again
155 }
156
157 if (doRestart)
158 {
159 restart();
160 }
161
162 Log << type() << " " << name() << " write:" << nl
163 << " Calculating averages" << nl;
164
165 forAll(faItems_, fieldi)
166 {
167 faItems_[fieldi].evolve(obr());
168 }
169
170 storeWindowFields<scalar>();
171 storeWindowFields<vector>();
172 storeWindowFields<sphericalTensor>();
173 storeWindowFields<symmTensor>();
174 storeWindowFields<tensor>();
175
176 addMeanSqrToPrime2Mean<scalar, scalar>();
177 addMeanSqrToPrime2Mean<vector, symmTensor>();
178
179 calculateMeanFields<scalar>();
180 calculateMeanFields<vector>();
181 calculateMeanFields<sphericalTensor>();
182 calculateMeanFields<symmTensor>();
183 calculateMeanFields<tensor>();
184
185 calculatePrime2MeanFields<scalar, scalar>();
186 calculatePrime2MeanFields<vector, symmTensor>();
187
188 Log << endl;
189}
190
191
193{
194 Log << " Writing average fields" << endl;
195
196 writeFields<scalar>();
197 writeFields<vector>();
198 writeFields<sphericalTensor>();
199 writeFields<symmTensor>();
200 writeFields<tensor>();
201
202 Log << endl;
203}
204
205
207{
208 for (const fieldAverageItem& item : faItems_)
209 {
211 item.writeState(propsDict);
212 setProperty(item.fieldName(), propsDict);
213 }
214}
215
216
218{
219 if (restartOnRestart_ || restartOnOutput_)
220 {
221 Info<< " Starting averaging at time "
222 << obr().time().timeOutputValue()
223 << nl;
224 }
225 else
226 {
227 Info<< " Restarting averaging for fields:" << nl;
228
229
230 for (fieldAverageItem& item : faItems_)
231 {
232 const word& fieldName = item.fieldName();
233 if (foundProperty(fieldName))
234 {
235 dictionary fieldDict;
236 getDict(fieldName, fieldDict);
237 item.readState(fieldDict);
238
239 if (item.allowRestart())
240 {
241 scalar userTotalTime =
242 obr().time().timeToUserTime(item.totalTime());
243
244 Info<< " " << fieldName
245 << ": iters = " << item.totalIter()
246 << " time = " << userTotalTime << nl;
247 }
248 else
249 {
250 item.clear(obr(), true);
251
252 Info<< " " << fieldName
253 << ": starting averaging at time "
254 << obr().time().timeOutputValue() << endl;
255 }
256 }
257 else
258 {
259 Info<< " " << fieldName
260 << ": starting averaging at time "
261 << obr().time().timeOutputValue() << endl;
262 }
263 }
264 }
265}
266
267
268// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
269
271(
272 const word& name,
273 const Time& runTime,
274 const dictionary& dict
275)
276:
278 prevTimeIndex_(-1),
279 initialised_(false),
280 restartOnRestart_(false),
281 restartOnOutput_(false),
282 periodicRestart_(false),
283 restartPeriod_(GREAT),
284 restartTime_(GREAT),
285 faItems_(),
286 periodIndex_(1)
287{
288 read(dict);
289}
290
291
292// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
293
295{
297
298 // Make certain that the values are consistent with the defaults:
299 initialised_ = false;
300 restartOnRestart_ = false;
301 restartOnOutput_ = false;
302 periodicRestart_ = false;
303 restartPeriod_ = GREAT;
304 restartTime_ = GREAT;
305
306 Info<< type() << " " << name() << ":" << nl;
307
308 dict.readIfPresent("restartOnRestart", restartOnRestart_);
309 dict.readIfPresent("restartOnOutput", restartOnOutput_);
310 dict.readIfPresent("periodicRestart", periodicRestart_);
311 dict.readEntry("fields", faItems_);
312
313 for (auto& item : faItems_)
314 {
315 item.setMeanFieldName(scopedName(item.meanFieldName()));
316 item.setPrime2MeanFieldName(scopedName(item.prime2MeanFieldName()));
317 }
318
319 const scalar currentTime = obr().time().value();
320
321 if (periodicRestart_)
322 {
323 scalar userRestartPeriod = dict.get<scalar>("restartPeriod");
324 restartPeriod_ = obr().time().userTimeToTime(userRestartPeriod);
325
326 if (restartPeriod_ > 0)
327 {
328 // Determine the appropriate interval for the next restart
329 periodIndex_ = 1;
330 while (currentTime > restartPeriod_*periodIndex_)
331 {
332 ++periodIndex_;
333 }
334
335 Info<< " Restart period " << userRestartPeriod
336 << " - next restart at " << (userRestartPeriod*periodIndex_)
337 << nl << endl;
338 }
339 else
340 {
341 periodicRestart_ = false;
342
343 Info<< " Restart period " << userRestartPeriod
344 << " - ignored"
345 << nl << endl;
346 }
347 }
348
349 scalar userRestartTime = 0;
350 if (dict.readIfPresent("restartTime", userRestartTime))
351 {
352 restartTime_ = obr().time().userTimeToTime(userRestartTime);
353
354 if (currentTime > restartTime_)
355 {
356 // The restart time is already in the past - ignore
357 restartTime_ = GREAT;
358 }
359 else
360 {
361 Info<< " Restart scheduled at time " << userRestartTime
362 << nl << endl;
363 }
364 }
365
366 readAveragingProperties();
367
368 Info<< endl;
369
370 return true;
371}
372
373
375{
376 calcAverages();
377
378 return true;
379}
380
381
383{
384 writeAverages();
385 writeAveragingProperties();
386
387 if (restartOnOutput_)
388 {
389 restart();
390 }
391
392 return true;
393}
394
395
396// ************************************************************************* //
#define Log
Definition: PDRblock.C:35
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
virtual bool read()
Re-read model coefficients if they have changed.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Abstract base-class for Time/database function objects.
const word & name() const noexcept
Return the name of this functionObject.
virtual const word & type() const =0
Runtime type information.
Helper class to describe what form of averaging to apply. A set will be applied to each base field in...
Computes ensemble- and/or time-based field averages, with optional windowing, for a user-specified se...
Definition: fieldAverage.H:262
label prevTimeIndex_
Time at last call, prevents repeated averaging.
Definition: fieldAverage.H:268
void writeAveragingProperties()
Write averaging properties - steps and time.
Definition: fieldAverage.C:206
void restart()
Restart averaging for restartOnOutput.
Definition: fieldAverage.C:104
void initialize()
Reset lists (clear existing values) and initialize averaging.
Definition: fieldAverage.C:48
virtual void calcAverages()
Main calculation routine.
Definition: fieldAverage.C:119
bool initialised_
Initialised flag.
Definition: fieldAverage.H:271
List< fieldAverageItem > faItems_
List of field average items, describing what averages to be.
Definition: fieldAverage.H:290
virtual void writeAverages() const
Write averages.
Definition: fieldAverage.C:192
void readAveragingProperties()
Read averaging properties - steps and time.
Definition: fieldAverage.C:217
virtual bool execute()
Calculate the field averages.
Definition: fieldAverage.C:374
virtual bool write()
Write the field averages.
Definition: fieldAverage.C:382
virtual bool read(const dictionary &)
Read the field average data.
Definition: fieldAverage.C:294
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
virtual const objectRegistry & obr() const
The region or sub-region registry being used.
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
engineTime & runTime
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for OpenFOAM.
messageStream Info
Information stream (stdout output on master, null elsewhere)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
IOdictionary propsDict(dictIO)