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-2015 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  faItems_[fieldI].active() = true;
37 
38  const word& fieldName = faItems_[fieldI].fieldName();
39  const word& meanFieldName = faItems_[fieldI].meanFieldName();
40 
41  Info<< " Reading/initialising field " << meanFieldName << endl;
42 
43  if (obr_.foundObject<Type>(meanFieldName))
44  {
45  // do nothing
46  }
47  else if (obr_.found(meanFieldName))
48  {
49  Info<< " Cannot allocate average field " << meanFieldName
50  << " since an object with that name already exists."
51  << " Disabling averaging for field." << endl;
52 
53  faItems_[fieldI].mean() = false;
54  }
55  else
56  {
57  const Type& baseField = obr_.lookupObject<Type>(fieldName);
58 
59  // Store on registry
60  obr_.store
61  (
62  new Type
63  (
64  IOobject
65  (
66  meanFieldName,
68  obr_,
73  ),
74  1*baseField
75  )
76  );
77  }
78 }
79 
80 
81 template<class Type>
83 {
84  if (faItems_[fieldI].mean())
85  {
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<surfFieldType>(fieldName))
96  {
97  addMeanFieldType<surfFieldType>(fieldI);
98  }
99  }
100 }
101 
102 
103 template<class Type1, class Type2>
105 {
106  const word& fieldName = faItems_[fieldI].fieldName();
107  const word& meanFieldName = faItems_[fieldI].meanFieldName();
108  const word& prime2MeanFieldName = faItems_[fieldI].prime2MeanFieldName();
109 
110  Info<< " Reading/initialising field " << prime2MeanFieldName << nl;
111 
112  if (obr_.foundObject<Type2>(prime2MeanFieldName))
113  {
114  // do nothing
115  }
116  else if (obr_.found(prime2MeanFieldName))
117  {
118  Info<< " 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;
156 
157  typedef GeometricField<Type2, fvPatchField, volMesh> volFieldType2;
159 
160  if (faItems_[fieldI].prime2Mean())
161  {
162  const word& fieldName = faItems_[fieldI].fieldName();
163 
164  if (!faItems_[fieldI].mean())
165  {
167  (
168  "void Foam::fieldAverage::addPrime2MeanField(const label) const"
169  )
170  << "To calculate the prime-squared average, the "
171  << "mean average must also be selected for field "
172  << fieldName << nl << exit(FatalError);
173  }
174 
175  if (obr_.foundObject<volFieldType1>(fieldName))
176  {
177  addPrime2MeanFieldType<volFieldType1, volFieldType2>(fieldI);
178  }
179  else if (obr_.foundObject<surfFieldType1>(fieldName))
180  {
181  addPrime2MeanFieldType<surfFieldType1, surfFieldType2>(fieldI);
182  }
183  }
184 }
185 
186 
187 template<class Type>
189 {
190  const word& fieldName = faItems_[fieldI].fieldName();
191 
192  if (obr_.foundObject<Type>(fieldName))
193  {
194  const Type& baseField = obr_.lookupObject<Type>(fieldName);
195 
196  Type& meanField = const_cast<Type&>
197  (
198  obr_.lookupObject<Type>(faItems_[fieldI].meanFieldName())
199  );
200 
201  scalar dt = obr_.time().deltaTValue();
202  scalar Dt = totalTime_[fieldI];
203 
204  if (faItems_[fieldI].iterBase())
205  {
206  dt = 1.0;
207  Dt = scalar(totalIter_[fieldI]);
208  }
209 
210  scalar alpha = (Dt - dt)/Dt;
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  alpha = (w - dt)/w;
220  beta = dt/w;
221  }
222  }
223 
224  meanField = alpha*meanField + beta*baseField;
225  }
226 }
227 
228 
229 template<class Type>
231 {
232  typedef GeometricField<Type, fvPatchField, volMesh> volFieldType;
234 
235  forAll(faItems_, i)
236  {
237  if (faItems_[i].mean())
238  {
239  calculateMeanFieldType<volFieldType>(i);
240  calculateMeanFieldType<surfFieldType>(i);
241  }
242  }
243 }
244 
245 
246 template<class Type1, class Type2>
248 {
249  const word& fieldName = faItems_[fieldI].fieldName();
250 
251  if (obr_.foundObject<Type1>(fieldName))
252  {
253  const Type1& baseField = obr_.lookupObject<Type1>(fieldName);
254  const Type1& meanField =
255  obr_.lookupObject<Type1>(faItems_[fieldI].meanFieldName());
256 
257  Type2& prime2MeanField = const_cast<Type2&>
258  (
259  obr_.lookupObject<Type2>(faItems_[fieldI].prime2MeanFieldName())
260  );
261 
262  scalar dt = obr_.time().deltaTValue();
263  scalar Dt = totalTime_[fieldI];
264 
265  if (faItems_[fieldI].iterBase())
266  {
267  dt = 1.0;
268  Dt = scalar(totalIter_[fieldI]);
269  }
270 
271  scalar alpha = (Dt - dt)/Dt;
272  scalar beta = dt/Dt;
273 
274  if (faItems_[fieldI].window() > 0)
275  {
276  const scalar w = faItems_[fieldI].window();
277 
278  if (Dt - dt >= w)
279  {
280  alpha = (w - dt)/w;
281  beta = dt/w;
282  }
283  }
284 
285  prime2MeanField =
286  alpha*prime2MeanField
287  + beta*sqr(baseField)
288  - sqr(meanField);
289  }
290 }
291 
292 
293 template<class Type1, class Type2>
295 {
296  typedef GeometricField<Type1, fvPatchField, volMesh> volFieldType1;
298 
299  typedef GeometricField<Type2, fvPatchField, volMesh> volFieldType2;
301 
302  forAll(faItems_, i)
303  {
304  if (faItems_[i].prime2Mean())
305  {
306  calculatePrime2MeanFieldType<volFieldType1, volFieldType2>(i);
307  calculatePrime2MeanFieldType<surfFieldType1, surfFieldType2>(i);
308  }
309  }
310 }
311 
312 
313 template<class Type1, class Type2>
315 {
316  const word& fieldName = faItems_[fieldI].fieldName();
317 
318  if (obr_.foundObject<Type1>(fieldName))
319  {
320  const Type1& meanField =
321  obr_.lookupObject<Type1>(faItems_[fieldI].meanFieldName());
322 
323  Type2& prime2MeanField = const_cast<Type2&>
324  (
325  obr_.lookupObject<Type2>(faItems_[fieldI].prime2MeanFieldName())
326  );
327 
328  prime2MeanField += sqr(meanField);
329  }
330 }
331 
332 
333 template<class Type1, class Type2>
335 {
336  typedef GeometricField<Type1, fvPatchField, volMesh> volFieldType1;
338 
339  typedef GeometricField<Type2, fvPatchField, volMesh> volFieldType2;
341 
342  forAll(faItems_, i)
343  {
344  if (faItems_[i].prime2Mean())
345  {
346  addMeanSqrToPrime2MeanType<volFieldType1, volFieldType2>(i);
347  addMeanSqrToPrime2MeanType<surfFieldType1, surfFieldType2>(i);
348  }
349  }
350 }
351 
352 
353 template<class Type>
354 void Foam::fieldAverage::writeFieldType(const word& fieldName) const
355 {
356  if (obr_.foundObject<Type>(fieldName))
357  {
358  const Type& f = obr_.lookupObject<Type>(fieldName);
359  f.write();
360  }
361 }
362 
363 
364 template<class Type>
366 {
367  typedef GeometricField<Type, fvPatchField, volMesh> volFieldType;
369 
370  forAll(faItems_, i)
371  {
372  if (faItems_[i].mean())
373  {
374  const word& fieldName = faItems_[i].meanFieldName();
375  writeFieldType<volFieldType>(fieldName);
376  writeFieldType<surfFieldType>(fieldName);
377  }
378  if (faItems_[i].prime2Mean())
379  {
380  const word& fieldName = faItems_[i].prime2MeanFieldName();
381  writeFieldType<volFieldType>(fieldName);
382  writeFieldType<surfFieldType>(fieldName);
383  }
384  }
385 }
386 
387 
388 // ************************************************************************* //
void addMeanSqrToPrime2Mean() const
Add mean-squared field value to prime-squared mean field.
void writeFields() const
Write fields.
const objectRegistry & obr_
Database this class is registered to.
Definition: fieldAverage.H:175
void calculatePrime2MeanFieldType(const label fieldI) const
Calculate prime-squared average fields.
Foam::surfaceFields.
List< label > totalIter_
Iteration steps counter.
Definition: fieldAverage.H:199
labelList f(nPoints)
void calculateMeanFields() const
Calculate mean average fields.
void calculatePrime2MeanFields() const
Calculate prime-squared average fields.
bool foundObject(const word &name) const
Is the named Type found?
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
void addPrime2MeanField(const label fieldI)
Add prime-squared average field to database.
void addPrime2MeanFieldType(const label fieldI)
Add prime-squared average field to database.
A class for handling words, derived from string.
Definition: word.H:59
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
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:741
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
void addMeanSqrToPrime2MeanType(const label fieldI) const
Add mean-squared field value to prime-squared mean field.
messageStream Info
virtual dimensionedScalar startTime() const
Return start time.
Definition: Time.C:872
List< fieldAverageItem > faItems_
List of field average items, describing what averages to be.
Definition: fieldAverage.H:194
static const char nl
Definition: Ostream.H:260
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
void writeFieldType(const word &fieldName) const
Write fields.
scalar deltaTValue() const
Return time step value.
Definition: TimeState.H:100
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
void store()
Transfer ownership of this object to its registry.
Definition: regIOobjectI.H:34
#define forAll(list, i)
Definition: UList.H:421
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
Switch resetOnOutput_
Reset the averaging process on output flag.
Definition: fieldAverage.H:187
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:314
Generic GeometricField class.
error FatalError
void calculateMeanFieldType(const label fieldI) const
Calculate mean average fields.
void addMeanFieldType(const label fieldI)
Add mean average field to database.
const Time & time() const
Return time.
void addMeanField(const label fieldI)
Add mean average field to database.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:109
const Type & value() const
Return const reference to value.
List< scalar > totalTime_
Total time counter.
Definition: fieldAverage.H:202