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