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-2022 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 Op>
72 (
73  const scalarField& values,
74  Result<scalar>& result,
75  const Op& op
76 ) const
77 {
78  label i = 0;
79  forAll(values, j)
80  {
81  if (op(values[j], values[i]))
82  {
83  i = j;
84  }
85  }
86 
87  result.value = values[i];
88  result.celli = isNull(cellIDs()) ? i : cellIDs()[i];
89  result.proci = Pstream::parRun() ? Pstream::myProcNo() : -1;
90  result.cc = fieldValue::mesh_.C()[result.celli];
91 
92  reduce
93  (
94  result,
95  [&op](const Result<scalar>& a, const Result<scalar>& b)
96  {
97  return op(a.value, b.value) ? a : b;
98  }
99  );
100 }
101 
102 
103 template<class Type, class ResultType>
105 (
106  const Field<Type>& values,
107  const scalarField& weights,
108  const scalarField& V,
109  Result<ResultType>& result
110 ) const
111 {
112  return false;
113 }
114 
115 
116 template<class Type>
118 (
119  const Field<Type>& values,
120  const scalarField& weights,
121  const scalarField& V,
122  Result<Type>& result
123 ) const
124 {
125  return processValuesTypeType(values, weights, V, result);
126 }
127 
128 
129 template<class Type>
131 (
132  const Field<Type>& values,
133  const scalarField& weights,
134  const scalarField& V,
135  Result<scalar>& result
136 ) const
137 {
138  switch (operation_)
139  {
140  case operationType::minMag:
141  {
142  compareScalars(mag(values), result, lessOp<scalar>());
143  return true;
144  }
145  case operationType::maxMag:
146  {
147  compareScalars(mag(values), result, greaterOp<scalar>());
148  return true;
149  }
150  default:
151  {
152  return false;
153  }
154  }
155 }
156 
157 
158 template<class Type>
160 (
161  const Field<Type>& values,
162  const scalarField& weights,
163  const scalarField& V,
164  Result<Type>& result
165 ) const
166 {
167  switch (operation_)
168  {
169  case operationType::sum:
170  {
171  result.value = gSum(weights*values);
172  return true;
173  }
175  {
176  result.value = gSum(cmptMag(values));
177  return true;
178  }
180  {
181  result.value = gSum(weights*values)/max(gSum(weights), vSmall);
182  return true;
183  }
184  case operationType::volAverage:
185  {
186  result.value = gSum(weights*V*values)/max(gSum(weights*V), vSmall);
187  return true;
188  }
189  case operationType::volIntegrate:
190  {
191  result.value = gSum(weights*V*values);
192  return true;
193  }
194  case operationType::min:
195  {
196  result.value = gMin(values);
197  return true;
198  }
199  case operationType::max:
200  {
201  result.value = gMax(values);
202  return true;
203  }
204  case operationType::CoV:
205  {
206  Type meanValue = gSum(values*V)/this->V();
207 
208  const label nComp = pTraits<Type>::nComponents;
209 
210  for (direction d=0; d<nComp; ++d)
211  {
212  scalarField vals(values.component(d));
213  scalar mean = component(meanValue, d);
214  scalar& res = setComponent(result.value, d);
215 
216  res = sqrt(gSum(V*sqr(vals - mean))/this->V())/mean;
217  }
218 
219  return true;
220  }
221  case operationType::none:
222  {
223  return true;
224  }
225  default:
226  {
227  return false;
228  }
229  }
230 }
231 
232 
233 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
234 
235 template<class Type>
237 (
238  const word& fieldName,
239  const scalarField& weights,
240  const scalarField& V
241 )
242 {
243  const bool ok = validField<Type>(fieldName);
244 
245  if (ok)
246  {
247  // Get the values
248  Field<Type> values(getFieldValues<Type>(fieldName));
249 
250  // Write raw values if specified
251  if (writeFields_)
252  {
254  (
255  IOobject
256  (
257  fieldName + '_' + regionTypeNames_[regionType_]
258  + '-' + volRegion::regionName_,
259  obr_.time().timeName(),
260  obr_,
261  IOobject::NO_READ,
262  IOobject::NO_WRITE
263  ),
264  (weights*values).ref()
265  ).write();
266  }
267 
268  // Do the operation
269  if (operation_ != operationType::none)
270  {
271  // Apply scale factor
272  values *= scaleFactor_;
273 
274  bool ok = false;
275 
276  #define writeValuesFieldType(fieldType, none) \
277  ok = \
278  ok \
279  || writeValues<Type, fieldType> \
280  ( \
281  fieldName, \
282  values, \
283  weights, \
284  V \
285  );
287  #undef writeValuesFieldType
288 
289  if (!ok)
290  {
292  << "Operation " << operationTypeNames_[operation_]
293  << " not available for values of type "
295  << exit(FatalError);
296  }
297  }
298  }
299 
300  return ok;
301 }
302 
303 
304 template<class Type, class ResultType>
306 (
307  const word& fieldName,
308  const Field<Type>& values,
309  const scalarField& weights,
310  const scalarField& V
311 )
312 {
313  Result<ResultType> result({Zero, -1, -1, point::uniform(NaN)});
314 
315  if (processValues(values, weights, V, result))
316  {
317  // Add to result dictionary, over-writing any previous entry
318  resultDict_.add(fieldName, result.value, true);
319 
320  if (Pstream::master())
321  {
322  file() << tab << result.value;
323 
324  Log << " " << operationTypeNames_[operation_]
325  << "(" << volRegion::regionName_ << ") of " << fieldName
326  << " = " << result.value;
327 
328  if (result.celli != -1)
329  {
330  Log << " at location " << result.cc;
331  if (writeLocation_) file() << tab << result.cc;
332  }
333 
334  if (result.celli != -1)
335  {
336  Log << " in cell " << result.celli;
337  if (writeLocation_) file() << tab << result.celli;
338  }
339 
340  if (result.proci != -1)
341  {
342  Log << " on processor " << result.proci;
343  if (writeLocation_) file() << tab << result.proci;
344  }
345 
346  Log << endl;
347  }
348 
349  return true;
350  }
351 
352  return false;
353 }
354 
355 
356 template<class Type>
359 (
360  const Field<Type>& field
361 ) const
362 {
363  if (isNull(cellIDs()))
364  {
365  return field;
366  }
367  else
368  {
369  return tmp<Field<Type>>(new Field<Type>(field, cellIDs()));
370  }
371 }
372 
373 
374 // ************************************************************************* //
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
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
void compareScalars(const scalarField &values, Result< scalar > &result, const Op &op) const
Apply a comparison operation to the values, returning the limiting.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const char tab
Definition: Ostream.H:259
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
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:51
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)
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)
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.
dimensionSet cmptMag(const dimensionSet &)
Definition: dimensionSet.C:275
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
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:455
tmp< Field< Type > > getFieldValues(const word &fieldName) const
Insert field values into values list.
Type gMax(const FieldField< Field, Type > &f)
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:98
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)