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-2018 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  else if (obr_.found(meanFieldName))
44  {
45  Log << " Cannot allocate average field " << meanFieldName
46  << " since an object with that name already exists."
47  << " Disabling averaging for field." << endl;
48 
49  faItems_[fieldi].mean() = false;
50  }
51  else
52  {
53  const Type& baseField = obr_.lookupObject<Type>(fieldName);
54 
55  // Store on registry
56  obr_.store
57  (
58  new Type
59  (
60  IOobject
61  (
62  meanFieldName,
64  obr_,
69  ),
70  1*baseField
71  )
72  );
73  }
74 }
75 
76 
77 template<class Type>
79 {
80  if (faItems_[fieldi].mean())
81  {
83  VolFieldType;
84 
86  SurfaceFieldType;
87 
88  const word& fieldName = faItems_[fieldi].fieldName();
89 
90  if (obr_.foundObject<VolFieldType>(fieldName))
91  {
92  addMeanFieldType<VolFieldType>(fieldi);
93  }
94  else if (obr_.foundObject<SurfaceFieldType>(fieldName))
95  {
96  addMeanFieldType<SurfaceFieldType>(fieldi);
97  }
98  }
99 }
100 
101 
102 template<class Type1, class Type2>
104 (
105  const label fieldi
106 )
107 {
108  const word& fieldName = faItems_[fieldi].fieldName();
109  const word& meanFieldName = faItems_[fieldi].meanFieldName();
110  const word& prime2MeanFieldName = faItems_[fieldi].prime2MeanFieldName();
111 
112  Log << " Reading/initialising field " << prime2MeanFieldName << nl;
113 
114  if (obr_.foundObject<Type2>(prime2MeanFieldName))
115  {}
116  else if (obr_.found(prime2MeanFieldName))
117  {
118  Log << " Cannot allocate average field " << prime2MeanFieldName
119  << " since an object with that name already exists."
120  << " Disabling averaging for field." << nl;
121 
122  faItems_[fieldi].prime2Mean() = false;
123  }
124  else
125  {
126  const Type1& baseField = obr_.lookupObject<Type1>(fieldName);
127  const Type1& meanField = obr_.lookupObject<Type1>(meanFieldName);
128 
129  // Store on registry
130  obr_.store
131  (
132  new Type2
133  (
134  IOobject
135  (
136  prime2MeanFieldName,
138  obr_,
143  ),
144  sqr(baseField) - sqr(meanField)
145  )
146  );
147  }
148 }
149 
150 
151 template<class Type1, class Type2>
153 {
154  typedef GeometricField<Type1, fvPatchField, volMesh> VolFieldType1;
155  typedef GeometricField<Type1, fvsPatchField, surfaceMesh> SurfaceFieldType1;
156 
157  typedef GeometricField<Type2, fvPatchField, volMesh> VolFieldType2;
158  typedef GeometricField<Type2, fvsPatchField, surfaceMesh> SurfaceFieldType2;
159 
160  if (faItems_[fieldi].prime2Mean())
161  {
162  const word& fieldName = faItems_[fieldi].fieldName();
163 
164  if (!faItems_[fieldi].mean())
165  {
167  << "To calculate the prime-squared average, the "
168  << "mean average must also be selected for field "
169  << fieldName << nl << exit(FatalError);
170  }
171 
172  if (obr_.foundObject<VolFieldType1>(fieldName))
173  {
174  addPrime2MeanFieldType<VolFieldType1, VolFieldType2>(fieldi);
175  }
176  else if (obr_.foundObject<SurfaceFieldType1>(fieldName))
177  {
178  addPrime2MeanFieldType<SurfaceFieldType1, SurfaceFieldType2>
179  (
180  fieldi
181  );
182  }
183  }
184 }
185 
186 
187 template<class Type>
189 (
190  const label fieldi
191 ) const
192 {
193  const word& fieldName = faItems_[fieldi].fieldName();
194 
195  if (obr_.foundObject<Type>(fieldName))
196  {
197  const Type& baseField = obr_.lookupObject<Type>(fieldName);
198 
199  Type& meanField =
200  obr_.lookupObjectRef<Type>(faItems_[fieldi].meanFieldName());
201 
202  scalar dt = obr_.time().deltaTValue();
203  scalar Dt = totalTime_[fieldi];
204 
205  if (faItems_[fieldi].iterBase())
206  {
207  dt = 1;
208  Dt = scalar(totalIter_[fieldi]);
209  }
210 
211  scalar beta = dt/Dt;
212 
213  if (faItems_[fieldi].window() > 0)
214  {
215  const scalar w = faItems_[fieldi].window();
216 
217  if (Dt - dt >= w)
218  {
219  beta = dt/w;
220  }
221  }
222 
223  meanField = (1 - beta)*meanField + beta*baseField;
224  }
225 }
226 
227 
228 template<class Type>
230 {
231  typedef GeometricField<Type, fvPatchField, volMesh> VolFieldType;
232  typedef GeometricField<Type, fvsPatchField, surfaceMesh> SurfaceFieldType;
233 
234  forAll(faItems_, i)
235  {
236  if (faItems_[i].mean())
237  {
238  calculateMeanFieldType<VolFieldType>(i);
239  calculateMeanFieldType<SurfaceFieldType>(i);
240  }
241  }
242 }
243 
244 
245 template<class Type1, class Type2>
247 (
248  const label fieldi
249 ) const
250 {
251  const word& fieldName = faItems_[fieldi].fieldName();
252 
253  if (obr_.foundObject<Type1>(fieldName))
254  {
255  const Type1& baseField = obr_.lookupObject<Type1>(fieldName);
256  const Type1& meanField =
257  obr_.lookupObject<Type1>(faItems_[fieldi].meanFieldName());
258 
259  Type2& prime2MeanField =
260  obr_.lookupObjectRef<Type2>(faItems_[fieldi].prime2MeanFieldName());
261 
262  scalar dt = obr_.time().deltaTValue();
263  scalar Dt = totalTime_[fieldi];
264 
265  if (faItems_[fieldi].iterBase())
266  {
267  dt = 1;
268  Dt = scalar(totalIter_[fieldi]);
269  }
270 
271  scalar beta = dt/Dt;
272 
273  if (faItems_[fieldi].window() > 0)
274  {
275  const scalar w = faItems_[fieldi].window();
276 
277  if (Dt - dt >= w)
278  {
279  beta = dt/w;
280  }
281  }
282 
283  prime2MeanField =
284  (1 - beta)*prime2MeanField
285  + beta*sqr(baseField)
286  - sqr(meanField);
287  }
288 }
289 
290 
291 template<class Type1, class Type2>
293 {
294  typedef GeometricField<Type1, fvPatchField, volMesh> VolFieldType1;
295  typedef GeometricField<Type1, fvsPatchField, surfaceMesh> SurfaceFieldType1;
296 
297  typedef GeometricField<Type2, fvPatchField, volMesh> VolFieldType2;
298  typedef GeometricField<Type2, fvsPatchField, surfaceMesh> SurfaceFieldType2;
299 
300  forAll(faItems_, i)
301  {
302  if (faItems_[i].prime2Mean())
303  {
304  calculatePrime2MeanFieldType<VolFieldType1, VolFieldType2>(i);
305  calculatePrime2MeanFieldType<SurfaceFieldType1, SurfaceFieldType2>
306  (
307  i
308  );
309  }
310  }
311 }
312 
313 
314 template<class Type1, class Type2>
316 (
317  const label fieldi
318 ) const
319 {
320  const word& fieldName = faItems_[fieldi].fieldName();
321 
322  if (obr_.foundObject<Type1>(fieldName))
323  {
324  const Type1& meanField =
325  obr_.lookupObject<Type1>(faItems_[fieldi].meanFieldName());
326 
327  Type2& prime2MeanField =
328  obr_.lookupObjectRef<Type2>(faItems_[fieldi].prime2MeanFieldName());
329 
330  prime2MeanField += sqr(meanField);
331  }
332 }
333 
334 
335 template<class Type1, class Type2>
337 {
338  typedef GeometricField<Type1, fvPatchField, volMesh> VolFieldType1;
339  typedef GeometricField<Type1, fvsPatchField, surfaceMesh> SurfaceFieldType1;
340 
341  typedef GeometricField<Type2, fvPatchField, volMesh> VolFieldType2;
342  typedef GeometricField<Type2, fvsPatchField, surfaceMesh> SurfaceFieldType2;
343 
344  forAll(faItems_, i)
345  {
346  if (faItems_[i].prime2Mean())
347  {
348  addMeanSqrToPrime2MeanType<VolFieldType1, VolFieldType2>(i);
349  addMeanSqrToPrime2MeanType<SurfaceFieldType1, SurfaceFieldType2>(i);
350  }
351  }
352 }
353 
354 
355 template<class Type>
357 (
358  const word& fieldName
359 ) const
360 {
361  if (obr_.foundObject<Type>(fieldName))
362  {
363  const Type& f = obr_.lookupObject<Type>(fieldName);
364  f.write();
365  }
366 }
367 
368 
369 template<class Type>
371 {
372  typedef GeometricField<Type, fvPatchField, volMesh> VolFieldType;
373  typedef GeometricField<Type, fvsPatchField, surfaceMesh> SurfaceFieldType;
374 
375  forAll(faItems_, i)
376  {
377  if (faItems_[i].mean())
378  {
379  const word& fieldName = faItems_[i].meanFieldName();
380  writeFieldType<VolFieldType>(fieldName);
381  writeFieldType<SurfaceFieldType>(fieldName);
382  }
383  if (faItems_[i].prime2Mean())
384  {
385  const word& fieldName = faItems_[i].prime2MeanFieldName();
386  writeFieldType<VolFieldType>(fieldName);
387  writeFieldType<SurfaceFieldType>(fieldName);
388  }
389  }
390 }
391 
392 
393 // ************************************************************************* //
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: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
error FatalError
#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:206
virtual dimensionedScalar startTime() const
Return start time.
Definition: Time.C:781
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:256
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:626
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:185
List< label > totalIter_
Iteration steps counter.
Definition: fieldAverage.H:203
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
List< fieldAverageItem > faItems_
List of field average items, describing what averages to be.
Definition: fieldAverage.H:198
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:265
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.
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
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