fieldAverageTemplates.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 "fieldAverageItem.H"
27 #include "volFields.H"
28 #include "surfaceFields.H"
29 #include "OFstream.H"
30 
31 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
32 
33 template<class Type>
35 {
36  const word& fieldName = faItems_[fieldi].fieldName();
37  const word& meanFieldName = faItems_[fieldi].meanFieldName();
38 
39  Log << " Reading/initialising field " << meanFieldName << endl;
40 
41  if (obr_.foundObject<Type>(meanFieldName))
42  {
43  // do nothing
44  }
45  else if (obr_.found(meanFieldName))
46  {
47  Log << " Cannot allocate average field " << meanFieldName
48  << " since an object with that name already exists."
49  << " Disabling averaging for field." << endl;
50 
51  faItems_[fieldi].mean() = false;
52  }
53  else
54  {
55  const Type& baseField = obr_.lookupObject<Type>(fieldName);
56 
57  // Store on registry
58  obr_.store
59  (
60  new Type
61  (
62  IOobject
63  (
64  meanFieldName,
66  obr_,
71  ),
72  1*baseField
73  )
74  );
75  }
76 }
77 
78 
79 template<class Type>
81 {
82  if (faItems_[fieldi].mean())
83  {
85  VolFieldType;
86 
88  SurfaceFieldType;
89 
90  const word& fieldName = faItems_[fieldi].fieldName();
91 
92  if (obr_.foundObject<VolFieldType>(fieldName))
93  {
94  addMeanFieldType<VolFieldType>(fieldi);
95  }
96  else if (obr_.foundObject<SurfaceFieldType>(fieldName))
97  {
98  addMeanFieldType<SurfaceFieldType>(fieldi);
99  }
100  }
101 }
102 
103 
104 template<class Type1, class Type2>
106 (
107  const label fieldi
108 )
109 {
110  const word& fieldName = faItems_[fieldi].fieldName();
111  const word& meanFieldName = faItems_[fieldi].meanFieldName();
112  const word& prime2MeanFieldName = faItems_[fieldi].prime2MeanFieldName();
113 
114  Log << " Reading/initialising field " << prime2MeanFieldName << nl;
115 
116  if (obr_.foundObject<Type2>(prime2MeanFieldName))
117  {
118  // do nothing
119  }
120  else if (obr_.found(prime2MeanFieldName))
121  {
122  Log << " Cannot allocate average field " << prime2MeanFieldName
123  << " since an object with that name already exists."
124  << " Disabling averaging for field." << nl;
125 
126  faItems_[fieldi].prime2Mean() = false;
127  }
128  else
129  {
130  const Type1& baseField = obr_.lookupObject<Type1>(fieldName);
131  const Type1& meanField = obr_.lookupObject<Type1>(meanFieldName);
132 
133  // Store on registry
134  obr_.store
135  (
136  new Type2
137  (
138  IOobject
139  (
140  prime2MeanFieldName,
142  obr_,
147  ),
148  sqr(baseField) - sqr(meanField)
149  )
150  );
151  }
152 }
153 
154 
155 template<class Type1, class Type2>
157 {
158  typedef GeometricField<Type1, fvPatchField, volMesh> VolFieldType1;
159  typedef GeometricField<Type1, fvsPatchField, surfaceMesh> SurfaceFieldType1;
160 
161  typedef GeometricField<Type2, fvPatchField, volMesh> VolFieldType2;
162  typedef GeometricField<Type2, fvsPatchField, surfaceMesh> SurfaceFieldType2;
163 
164  if (faItems_[fieldi].prime2Mean())
165  {
166  const word& fieldName = faItems_[fieldi].fieldName();
167 
168  if (!faItems_[fieldi].mean())
169  {
171  << "To calculate the prime-squared average, the "
172  << "mean average must also be selected for field "
173  << fieldName << nl << exit(FatalError);
174  }
175 
176  if (obr_.foundObject<VolFieldType1>(fieldName))
177  {
178  addPrime2MeanFieldType<VolFieldType1, VolFieldType2>(fieldi);
179  }
180  else if (obr_.foundObject<SurfaceFieldType1>(fieldName))
181  {
182  addPrime2MeanFieldType<SurfaceFieldType1, SurfaceFieldType2>
183  (
184  fieldi
185  );
186  }
187  }
188 }
189 
190 
191 template<class Type>
193 (
194  const label fieldi
195 ) const
196 {
197  const word& fieldName = faItems_[fieldi].fieldName();
198 
199  if (obr_.foundObject<Type>(fieldName))
200  {
201  const Type& baseField = obr_.lookupObject<Type>(fieldName);
202 
203  Type& meanField = const_cast<Type&>
204  (
205  obr_.lookupObject<Type>(faItems_[fieldi].meanFieldName())
206  );
207 
208  scalar dt = obr_.time().deltaTValue();
209  scalar Dt = totalTime_[fieldi];
210 
211  if (faItems_[fieldi].iterBase())
212  {
213  dt = 1.0;
214  Dt = scalar(totalIter_[fieldi]);
215  }
216 
217  scalar alpha = (Dt - dt)/Dt;
218  scalar beta = dt/Dt;
219 
220  if (faItems_[fieldi].window() > 0)
221  {
222  const scalar w = faItems_[fieldi].window();
223 
224  if (Dt - dt >= w)
225  {
226  alpha = (w - dt)/w;
227  beta = dt/w;
228  }
229  }
230 
231  meanField = alpha*meanField + beta*baseField;
232  }
233 }
234 
235 
236 template<class Type>
238 {
239  typedef GeometricField<Type, fvPatchField, volMesh> VolFieldType;
240  typedef GeometricField<Type, fvsPatchField, surfaceMesh> SurfaceFieldType;
241 
242  forAll(faItems_, i)
243  {
244  if (faItems_[i].mean())
245  {
246  calculateMeanFieldType<VolFieldType>(i);
247  calculateMeanFieldType<SurfaceFieldType>(i);
248  }
249  }
250 }
251 
252 
253 template<class Type1, class Type2>
255 (
256  const label fieldi
257 ) const
258 {
259  const word& fieldName = faItems_[fieldi].fieldName();
260 
261  if (obr_.foundObject<Type1>(fieldName))
262  {
263  const Type1& baseField = obr_.lookupObject<Type1>(fieldName);
264  const Type1& meanField =
265  obr_.lookupObject<Type1>(faItems_[fieldi].meanFieldName());
266 
267  Type2& prime2MeanField = const_cast<Type2&>
268  (
269  obr_.lookupObject<Type2>(faItems_[fieldi].prime2MeanFieldName())
270  );
271 
272  scalar dt = obr_.time().deltaTValue();
273  scalar Dt = totalTime_[fieldi];
274 
275  if (faItems_[fieldi].iterBase())
276  {
277  dt = 1.0;
278  Dt = scalar(totalIter_[fieldi]);
279  }
280 
281  scalar alpha = (Dt - dt)/Dt;
282  scalar beta = dt/Dt;
283 
284  if (faItems_[fieldi].window() > 0)
285  {
286  const scalar w = faItems_[fieldi].window();
287 
288  if (Dt - dt >= w)
289  {
290  alpha = (w - dt)/w;
291  beta = dt/w;
292  }
293  }
294 
295  prime2MeanField =
296  alpha*prime2MeanField
297  + beta*sqr(baseField)
298  - sqr(meanField);
299  }
300 }
301 
302 
303 template<class Type1, class Type2>
305 {
306  typedef GeometricField<Type1, fvPatchField, volMesh> VolFieldType1;
307  typedef GeometricField<Type1, fvsPatchField, surfaceMesh> SurfaceFieldType1;
308 
309  typedef GeometricField<Type2, fvPatchField, volMesh> VolFieldType2;
310  typedef GeometricField<Type2, fvsPatchField, surfaceMesh> SurfaceFieldType2;
311 
312  forAll(faItems_, i)
313  {
314  if (faItems_[i].prime2Mean())
315  {
316  calculatePrime2MeanFieldType<VolFieldType1, VolFieldType2>
317  (
318  i
319  );
320 
321  calculatePrime2MeanFieldType<SurfaceFieldType1, SurfaceFieldType2>
322  (
323  i
324  );
325  }
326  }
327 }
328 
329 
330 template<class Type1, class Type2>
332 (
333  const label fieldi
334 ) const
335 {
336  const word& fieldName = faItems_[fieldi].fieldName();
337 
338  if (obr_.foundObject<Type1>(fieldName))
339  {
340  const Type1& meanField =
341  obr_.lookupObject<Type1>(faItems_[fieldi].meanFieldName());
342 
343  Type2& prime2MeanField = const_cast<Type2&>
344  (
345  obr_.lookupObject<Type2>(faItems_[fieldi].prime2MeanFieldName())
346  );
347 
348  prime2MeanField += sqr(meanField);
349  }
350 }
351 
352 
353 template<class Type1, class Type2>
355 {
356  typedef GeometricField<Type1, fvPatchField, volMesh> VolFieldType1;
357  typedef GeometricField<Type1, fvsPatchField, surfaceMesh> SurfaceFieldType1;
358 
359  typedef GeometricField<Type2, fvPatchField, volMesh> VolFieldType2;
360  typedef GeometricField<Type2, fvsPatchField, surfaceMesh> SurfaceFieldType2;
361 
362  forAll(faItems_, i)
363  {
364  if (faItems_[i].prime2Mean())
365  {
366  addMeanSqrToPrime2MeanType<VolFieldType1, VolFieldType2>(i);
367  addMeanSqrToPrime2MeanType<SurfaceFieldType1, SurfaceFieldType2>(i);
368  }
369  }
370 }
371 
372 
373 template<class Type>
375 (
376  const word& fieldName
377 ) const
378 {
379  if (obr_.foundObject<Type>(fieldName))
380  {
381  const Type& f = obr_.lookupObject<Type>(fieldName);
382  f.write();
383  }
384 }
385 
386 
387 template<class Type>
389 {
390  typedef GeometricField<Type, fvPatchField, volMesh> VolFieldType;
391  typedef GeometricField<Type, fvsPatchField, surfaceMesh> SurfaceFieldType;
392 
393  forAll(faItems_, i)
394  {
395  if (faItems_[i].mean())
396  {
397  const word& fieldName = faItems_[i].meanFieldName();
398  writeFieldType<VolFieldType>(fieldName);
399  writeFieldType<SurfaceFieldType>(fieldName);
400  }
401  if (faItems_[i].prime2Mean())
402  {
403  const word& fieldName = faItems_[i].prime2MeanFieldName();
404  writeFieldType<VolFieldType>(fieldName);
405  writeFieldType<SurfaceFieldType>(fieldName);
406  }
407  }
408 }
409 
410 
411 // ************************************************************************* //
Foam::surfaceFields.
void addPrime2MeanFieldType(const label fieldi)
Add prime-squared average field to database.
const Time & time() const
Return time.
void addMeanField(const label fieldi)
Add mean average field to database.
#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
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
void calculateMeanFields() const
Calculate mean average fields.
error FatalError
bool foundObject(const word &name) const
Is the named Type found?
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
dimensionedSymmTensor sqr(const dimensionedVector &dv)
List< scalar > totalTime_
Total time counter.
Definition: fieldAverage.H:211
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
virtual dimensionedScalar startTime() const
Return start time.
Definition: Time.C:835
void addMeanFieldType(const label fieldi)
Add mean average field to database.
const Type & value() const
Return const reference to value.
Generic GeometricField class.
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:715
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
Switch restartOnOutput_
Restart the averaging process on output.
Definition: fieldAverage.H:190
List< label > totalIter_
Iteration steps counter.
Definition: fieldAverage.H:208
void writeFields() const
Write fields.
A class for handling words, derived from string.
Definition: word.H:59
List< fieldAverageItem > faItems_
List of field average items, describing what averages to be.
Definition: fieldAverage.H:203
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:109
void store()
Transfer ownership of this object to its registry.
Definition: regIOobjectI.H:34
void calculatePrime2MeanFieldType(const label fieldi) const
Calculate prime-squared average fields.
static const char nl
Definition: Ostream.H:262
labelList f(nPoints)
scalar deltaTValue() const
Return time step value.
Definition: TimeStateI.H:41
void addMeanSqrToPrime2Mean() const
Add mean-squared field value to prime-squared mean field.
void addPrime2MeanField(const label fieldi)
Add prime-squared average field to database.
const objectRegistry & obr_
Reference to the region objectRegistry.
void writeFieldType(const word &fieldName) const
Write fields.
void addMeanSqrToPrime2MeanType(const label fieldi) const
Add mean-squared field value to prime-squared mean field.
#define Log
Report write to Foam::Info if the local log switch is true.
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
void calculateMeanFieldType(const label fieldi) const
Calculate mean average fields.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
void calculatePrime2MeanFields() const
Calculate prime-squared average fields.