volFieldValueTemplates.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 "volFieldValue.H"
27 #include "volFields.H"
28 
29 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
30 
31 template<class Type>
33 (
34  const word& fieldName
35 ) const
36 {
38 
39  if (obr_.foundObject<vf>(fieldName))
40  {
41  return true;
42  }
43 
44  return false;
45 }
46 
47 
48 template<class Type>
51 (
52  const word& fieldName
53 ) const
54 {
56 
57  if (obr_.foundObject<vf>(fieldName))
58  {
59  return filterField(obr_.lookupObject<vf>(fieldName));
60  }
61 
63  << "Field " << fieldName << " not found in database"
64  << abort(FatalError);
65 
66  return tmp<Field<Type>>(nullptr);
67 }
68 
69 
70 template<class Type, class Op>
72 (
73  const Field<Type>& values,
74  Result<scalar>& result,
75  const Op& op
76 ) const
77 {
78  const scalarField magValues(mag(values));
79 
80  label i = 0;
81  forAll(magValues, j)
82  {
83  if (op(magValues[j], magValues[i]))
84  {
85  i = j;
86  }
87  }
88 
89  result.value = magValues[i];
90  result.celli = isNull(cellIDs()) ? i : cellIDs()[i];
91  result.proci = Pstream::parRun() ? Pstream::myProcNo() : -1;
92  result.cc = fieldValue::mesh_.C()[result.celli];
93 
94  reduce
95  (
96  result,
97  [&op](const Result<scalar>& a, const Result<scalar>& b)
98  {
99  return op(a.value, b.value) ? a : b;
100  }
101  );
102 }
103 
104 
105 template<class Type, class ResultType>
107 (
108  const Field<Type>& values,
109  const scalarField& weights,
110  const scalarField& V,
111  Result<ResultType>& result
112 ) const
113 {
114  return false;
115 }
116 
117 
118 template<class Type>
120 (
121  const Field<Type>& values,
122  const scalarField& weights,
123  const scalarField& V,
124  Result<Type>& result
125 ) const
126 {
127  return processValuesTypeType(values, weights, V, result);
128 }
129 
130 
131 template<class Type>
133 (
134  const Field<Type>& values,
135  const scalarField& weights,
136  const scalarField& V,
137  Result<scalar>& result
138 ) const
139 {
140  switch (operation_)
141  {
142  case operationType::minMag:
143  {
144  opMag(values, result, lessOp<scalar>());
145  return true;
146  }
147  case operationType::maxMag:
148  {
149  opMag(values, result, greaterOp<scalar>());
150  return true;
151  }
152  default:
153  {
154  return false;
155  }
156  }
157 }
158 
159 
160 template<class Type>
162 (
163  const Field<Type>& values,
164  const scalarField& weights,
165  const scalarField& V,
166  Result<Type>& result
167 ) const
168 {
169  switch (operation_)
170  {
171  case operationType::sum:
172  {
173  result.value = gSum(weights*values);
174  return true;
175  }
177  {
178  result.value = gSum(cmptMag(values));
179  return true;
180  }
182  {
183  result.value = gSum(weights*values)/max(gSum(weights), vSmall);
184  return true;
185  }
186  case operationType::volAverage:
187  {
188  result.value = gSum(weights*V*values)/max(gSum(weights*V), vSmall);
189  return true;
190  }
191  case operationType::volIntegrate:
192  {
193  result.value = gSum(weights*V*values);
194  return true;
195  }
196  case operationType::min:
197  {
198  result.value = gMin(values);
199  return true;
200  }
201  case operationType::max:
202  {
203  result.value = gMax(values);
204  return true;
205  }
206  case operationType::CoV:
207  {
208  Type meanValue = gSum(values*V)/this->V();
209 
210  const label nComp = pTraits<Type>::nComponents;
211 
212  for (direction d=0; d<nComp; ++d)
213  {
214  scalarField vals(values.component(d));
215  scalar mean = component(meanValue, d);
216  scalar& res = setComponent(result.value, d);
217 
218  res = sqrt(gSum(V*sqr(vals - mean))/this->V())/mean;
219  }
220 
221  return true;
222  }
223  case operationType::none:
224  {
225  return true;
226  }
227  default:
228  {
229  return false;
230  }
231  }
232 }
233 
234 
235 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
236 
237 template<class Type>
239 (
240  const word& fieldName,
241  const scalarField& weights,
242  const scalarField& V
243 )
244 {
245  const bool ok = validField<Type>(fieldName);
246 
247  if (ok)
248  {
249  // Get the values
250  Field<Type> values(getFieldValues<Type>(fieldName));
251 
252  // Write raw values if specified
253  if (writeFields_)
254  {
256  (
257  IOobject
258  (
259  fieldName + '_' + regionTypeNames_[regionType_]
260  + '-' + volRegion::regionName_,
261  obr_.time().timeName(),
262  obr_,
263  IOobject::NO_READ,
264  IOobject::NO_WRITE
265  ),
266  (weights*values).ref()
267  ).write();
268  }
269 
270  // Do the operation
271  if (operation_ != operationType::none)
272  {
273  // Apply scale factor
274  values *= scaleFactor_;
275 
276  bool ok = false;
277 
278  #define writeValuesFieldType(fieldType, none) \
279  ok = \
280  ok \
281  || writeValues<Type, fieldType> \
282  ( \
283  fieldName, \
284  values, \
285  weights, \
286  V \
287  );
289  #undef writeValuesFieldType
290 
291  if (!ok)
292  {
294  << "Operation " << operationTypeNames_[operation_]
295  << " not available for values of type "
297  << exit(FatalError);
298  }
299  }
300  }
301 
302  return ok;
303 }
304 
305 
306 template<class Type, class ResultType>
308 (
309  const word& fieldName,
310  const Field<Type>& values,
311  const scalarField& weights,
312  const scalarField& V
313 )
314 {
315  Result<ResultType> result({Zero, -1, -1, point::uniform(NaN)});
316 
317  if (processValues(values, weights, V, result))
318  {
319  // Add to result dictionary, over-writing any previous entry
320  resultDict_.add(fieldName, result.value, true);
321 
322  if (Pstream::master())
323  {
324  file() << tab << result.value;
325 
326  Log << " " << operationTypeNames_[operation_]
327  << "(" << volRegion::regionName_ << ") of " << fieldName
328  << " = " << result.value;
329 
330  if (result.celli != -1)
331  {
332  Log << " at location " << result.cc;
333  if (writeLocation_) file() << tab << result.cc;
334  }
335 
336  if (result.celli != -1)
337  {
338  Log << " in cell " << result.celli;
339  if (writeLocation_) file() << tab << result.celli;
340  }
341 
342  if (result.proci != -1)
343  {
344  Log << " on processor " << result.proci;
345  if (writeLocation_) file() << tab << result.proci;
346  }
347 
348  Log << endl;
349  }
350 
351  return true;
352  }
353 
354  return false;
355 }
356 
357 
358 template<class Type>
361 (
362  const Field<Type>& field
363 ) const
364 {
365  if (isNull(cellIDs()))
366  {
367  return field;
368  }
369  else
370  {
371  return tmp<Field<Type>>(new Field<Type>(field, cellIDs()));
372  }
373 }
374 
375 
376 // ************************************************************************* //
bool processValues(const Field< Type > &values, const scalarField &weights, const scalarField &V, Result< ResultType > &result) const
Apply the operation to the values, and return true if successful.
#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
static const char tab
Definition: Ostream.H:259
error FatalError
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:323
Type gMin(const FieldField< Field, Type > &f)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
uint8_t direction
Definition: direction.H:45
dimensionedScalar sqrt(const dimensionedScalar &ds)
const dimensionedScalar b
Wien displacement law constant: default SI units: [m K].
Definition: createFields.H:27
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Traits class for primitives.
Definition: pTraits.H:50
#define Op(opName, op)
Definition: ops.H:101
#define FOR_ALL_FIELD_TYPES(Macro,...)
Definition: fieldTypes.H:50
Generic GeometricField class.
bool processValuesTypeType(const Field< Type > &values, const scalarField &weights, const scalarField &V, Result< Type > &result) const
Apply a Type -> Type operation to the values.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
void opMag(const Field< Type > &values, Result< scalar > &result, const Op &op) const
Apply a comparison operation (min/max) to the field magnitude,.
bool isNull(const T &t)
Return true if t is a reference to the nullObject of type T.
Definition: nullObjectI.H:46
Type gSum(const FieldField< Field, Type > &f)
Pre-declare SubField and related Field type.
Definition: Field.H:56
bool validField(const word &fieldName) const
Return true if the field name is valid.
A class for handling words, derived from string.
Definition: word.H:59
tmp< Field< Type > > filterField(const Field< Type > &field) const
Filter a field according to cellIds.
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
static const zero Zero
Definition: zero.H:97
bool writeValues(const word &fieldName, const scalarField &weights, const scalarField &V)
Templated helper function to output field values.
errorManip< error > abort(error &err)
Definition: errorManip.H:131
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
Definition: Field.C:441
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
tmp< Field< Type > > getFieldValues(const word &fieldName) const
Insert field values into values list.
Type gMax(const FieldField< Field, Type > &f)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
#define Log
Report write to Foam::Info if the local log switch is true.
tmp< GeometricField< Type, fvPatchField, volMesh > > average(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Area-weighted average a surfaceField creating a volField.
Definition: fvcAverage.C:46
dimensioned< scalar > mag(const dimensioned< Type > &)
rDeltaTY field()
A class for managing temporary objects.
Definition: PtrList.H:53
label & setComponent(label &l, const direction)
Definition: label.H:86
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
#define writeValuesFieldType(fieldType, none)
A primitive field of type <T> with automated input and output.
Definition: IOField.H:50
dimensioned< scalar > sumMag(const DimensionedField< Type, GeoMesh > &df)