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-2021 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 
84  typedef typename VolFieldType::Internal InternalType;
85 
87  SurfaceFieldType;
88 
89  const word& fieldName = faItems_[fieldi].fieldName();
90 
91  if (obr_.foundObject<VolFieldType>(fieldName))
92  {
93  addMeanFieldType<VolFieldType>(fieldi);
94  }
95  else if (obr_.foundObject<InternalType>(fieldName))
96  {
97  addMeanFieldType<InternalType>(fieldi);
98  }
99  else if (obr_.foundObject<SurfaceFieldType>(fieldName))
100  {
101  addMeanFieldType<SurfaceFieldType>(fieldi);
102  }
103  }
104 }
105 
106 
107 template<class Type1, class Type2>
109 (
110  const label fieldi
111 )
112 {
113  const word& fieldName = faItems_[fieldi].fieldName();
114  const word& meanFieldName = faItems_[fieldi].meanFieldName();
115  const word& prime2MeanFieldName = faItems_[fieldi].prime2MeanFieldName();
116 
117  Log << " Reading/initialising field " << prime2MeanFieldName << nl;
118 
119  if (obr_.foundObject<Type2>(prime2MeanFieldName))
120  {}
121  else if (obr_.found(prime2MeanFieldName))
122  {
123  Log << " Cannot allocate average field " << prime2MeanFieldName
124  << " since an object with that name already exists."
125  << " Disabling averaging for field." << nl;
126 
127  faItems_[fieldi].prime2Mean() = false;
128  }
129  else
130  {
131  const Type1& baseField = obr_.lookupObject<Type1>(fieldName);
132  const Type1& meanField = obr_.lookupObject<Type1>(meanFieldName);
133 
134  // Store on registry
135  obr_.store
136  (
137  new Type2
138  (
139  IOobject
140  (
141  prime2MeanFieldName,
143  obr_,
148  ),
149  sqr(baseField) - sqr(meanField)
150  )
151  );
152  }
153 }
154 
155 
156 template<class Type1, class Type2>
158 {
159  typedef GeometricField<Type1, fvPatchField, volMesh> VolFieldType1;
160  typedef typename VolFieldType1::Internal InternalType1;
161  typedef GeometricField<Type1, fvsPatchField, surfaceMesh> SurfaceFieldType1;
162 
163  typedef GeometricField<Type2, fvPatchField, volMesh> VolFieldType2;
164  typedef typename VolFieldType2::Internal InternalType2;
165  typedef GeometricField<Type2, fvsPatchField, surfaceMesh> SurfaceFieldType2;
166 
167  if (faItems_[fieldi].prime2Mean())
168  {
169  const word& fieldName = faItems_[fieldi].fieldName();
170 
171  if (!faItems_[fieldi].mean())
172  {
174  << "To calculate the prime-squared average, the "
175  << "mean average must also be selected for field "
176  << fieldName << nl << exit(FatalError);
177  }
178 
179  if (obr_.foundObject<VolFieldType1>(fieldName))
180  {
181  addPrime2MeanFieldType<VolFieldType1, VolFieldType2>(fieldi);
182  }
183  else if (obr_.foundObject<InternalType1>(fieldName))
184  {
185  addPrime2MeanFieldType<InternalType1, InternalType2>(fieldi);
186  }
187  else if (obr_.foundObject<SurfaceFieldType1>(fieldName))
188  {
189  addPrime2MeanFieldType<SurfaceFieldType1, SurfaceFieldType2>
190  (
191  fieldi
192  );
193  }
194  }
195 }
196 
197 
198 template<class Type>
200 (
201  const label fieldi
202 ) const
203 {
204  const word& fieldName = faItems_[fieldi].fieldName();
205 
206  if (obr_.foundObject<Type>(fieldName))
207  {
208  const Type& baseField = obr_.lookupObject<Type>(fieldName);
209 
210  Type& meanField =
211  obr_.lookupObjectRef<Type>(faItems_[fieldi].meanFieldName());
212 
213  scalar dt = obr_.time().deltaTValue();
214  scalar Dt = totalTime_[fieldi];
215 
216  if (iterBase())
217  {
218  dt = 1;
219  Dt = scalar(totalIter_[fieldi]);
220  }
221 
222  scalar beta = dt/Dt;
223 
224  if (window() > 0)
225  {
226  const scalar w = window();
227 
228  if (Dt - dt >= w)
229  {
230  beta = dt/w;
231  }
232  }
233 
234  meanField = (1 - beta)*meanField + beta*baseField;
235  }
236 }
237 
238 
239 template<class Type>
241 {
242  typedef GeometricField<Type, fvPatchField, volMesh> VolFieldType;
243  typedef typename VolFieldType::Internal InternalType;
244  typedef GeometricField<Type, fvsPatchField, surfaceMesh> SurfaceFieldType;
245 
246  forAll(faItems_, i)
247  {
248  if (faItems_[i].mean())
249  {
250  calculateMeanFieldType<VolFieldType>(i);
251  calculateMeanFieldType<InternalType>(i);
252  calculateMeanFieldType<SurfaceFieldType>(i);
253  }
254  }
255 }
256 
257 
258 template<class Type1, class Type2>
260 (
261  const label fieldi
262 ) const
263 {
264  const word& fieldName = faItems_[fieldi].fieldName();
265 
266  if (obr_.foundObject<Type1>(fieldName))
267  {
268  const Type1& baseField = obr_.lookupObject<Type1>(fieldName);
269  const Type1& meanField =
270  obr_.lookupObject<Type1>(faItems_[fieldi].meanFieldName());
271 
272  Type2& prime2MeanField =
273  obr_.lookupObjectRef<Type2>(faItems_[fieldi].prime2MeanFieldName());
274 
275  scalar dt = obr_.time().deltaTValue();
276  scalar Dt = totalTime_[fieldi];
277 
278  if (iterBase())
279  {
280  dt = 1;
281  Dt = scalar(totalIter_[fieldi]);
282  }
283 
284  scalar beta = dt/Dt;
285 
286  if (window() > 0)
287  {
288  const scalar w = window();
289 
290  if (Dt - dt >= w)
291  {
292  beta = dt/w;
293  }
294  }
295 
296  prime2MeanField =
297  (1 - beta)*prime2MeanField
298  + beta*sqr(baseField)
299  - sqr(meanField);
300  }
301 }
302 
303 
304 template<class Type1, class Type2>
306 {
307  typedef GeometricField<Type1, fvPatchField, volMesh> VolFieldType1;
308  typedef typename VolFieldType1::Internal InternalType1;
309  typedef GeometricField<Type1, fvsPatchField, surfaceMesh> SurfaceFieldType1;
310 
311  typedef GeometricField<Type2, fvPatchField, volMesh> VolFieldType2;
312  typedef typename VolFieldType2::Internal InternalType2;
313  typedef GeometricField<Type2, fvsPatchField, surfaceMesh> SurfaceFieldType2;
314 
315  forAll(faItems_, i)
316  {
317  if (faItems_[i].prime2Mean())
318  {
319  calculatePrime2MeanFieldType<VolFieldType1, VolFieldType2>(i);
320  calculatePrime2MeanFieldType<InternalType1, InternalType2>(i);
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 =
344  obr_.lookupObjectRef<Type2>(faItems_[fieldi].prime2MeanFieldName());
345 
346  prime2MeanField += sqr(meanField);
347  }
348 }
349 
350 
351 template<class Type1, class Type2>
353 {
354  typedef GeometricField<Type1, fvPatchField, volMesh> VolFieldType1;
355  typedef typename VolFieldType1::Internal InternalType1;
356  typedef GeometricField<Type1, fvsPatchField, surfaceMesh> SurfaceFieldType1;
357 
358  typedef GeometricField<Type2, fvPatchField, volMesh> VolFieldType2;
359  typedef typename VolFieldType2::Internal InternalType2;
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<InternalType1, InternalType2>(i);
368  addMeanSqrToPrime2MeanType<SurfaceFieldType1, SurfaceFieldType2>(i);
369  }
370  }
371 }
372 
373 
374 template<class Type>
376 (
377  const word& fieldName
378 ) const
379 {
380  if (obr_.foundObject<Type>(fieldName))
381  {
382  const Type& f = obr_.lookupObject<Type>(fieldName);
383  f.write();
384  }
385 }
386 
387 
388 template<class Type>
390 {
391  typedef GeometricField<Type, fvPatchField, volMesh> VolFieldType;
392  typedef typename VolFieldType::Internal InternalType;
393  typedef GeometricField<Type, fvsPatchField, surfaceMesh> SurfaceFieldType;
394 
395  forAll(faItems_, i)
396  {
397  if (faItems_[i].mean())
398  {
399  const word& fieldName = faItems_[i].meanFieldName();
400  writeFieldType<VolFieldType>(fieldName);
401  writeFieldType<InternalType>(fieldName);
402  writeFieldType<SurfaceFieldType>(fieldName);
403  }
404  if (faItems_[i].prime2Mean())
405  {
406  const word& fieldName = faItems_[i].prime2MeanFieldName();
407  writeFieldType<VolFieldType>(fieldName);
408  writeFieldType<InternalType>(fieldName);
409  writeFieldType<SurfaceFieldType>(fieldName);
410  }
411  }
412 }
413 
414 
415 // ************************************************************************* //
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:323
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:797
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:636
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