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-2023 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 {
37  if
38  (
39  obr_.foundObject<VolField<Type>>(fieldName)
41  )
42  {
43  return true;
44  }
45 
46  return false;
47 }
48 
49 
50 template<class Type>
53 (
54  const word& fieldName
55 ) const
56 {
57  if (obr_.foundObject<VolField<Type>>(fieldName))
58  {
59  return filterField(obr_.lookupObject<VolField<Type>>(fieldName));
60  }
61 
62  if (obr_.foundObject<VolInternalField<Type>>(fieldName))
63  {
64  return filterField
65  (
66  obr_.lookupObject<VolInternalField<Type>>(fieldName)
67  );
68  }
69 
71  << "Field " << fieldName << " not found in database"
72  << abort(FatalError);
73 
74  return tmp<Field<Type>>(nullptr);
75 }
76 
77 
78 template<class Op>
80 (
81  const scalarField& values,
82  Result<scalar>& result,
83  const Op& op
84 ) const
85 {
86  label i = 0;
87  forAll(values, j)
88  {
89  if (op(values[j], values[i]))
90  {
91  i = j;
92  }
93  }
94 
95  result.value = values[i];
96  result.celli = celli(i);
97  result.proci = Pstream::parRun() ? Pstream::myProcNo() : -1;
98  result.cc = fieldValue::mesh_.C()[result.celli];
99 
100  reduce
101  (
102  result,
103  [&op](const Result<scalar>& a, const Result<scalar>& b)
104  {
105  return op(a.value, b.value) ? a : b;
106  }
107  );
108 }
109 
110 
111 template<class Type, class ResultType>
113 (
114  const Field<Type>& values,
115  const scalarField& weights,
116  const scalarField& V,
117  Result<ResultType>& result
118 ) const
119 {
120  return false;
121 }
122 
123 
124 template<class Type>
126 (
127  const Field<Type>& values,
128  const scalarField& weights,
129  const scalarField& V,
130  Result<Type>& result
131 ) const
132 {
133  return processValuesTypeType(values, weights, V, result);
134 }
135 
136 
137 template<class Type>
139 (
140  const Field<Type>& values,
141  const scalarField& weights,
142  const scalarField& V,
143  Result<scalar>& result
144 ) const
145 {
146  switch (operation_)
147  {
148  case operationType::minMag:
149  {
150  compareScalars(mag(values), result, lessOp<scalar>());
151  return true;
152  }
153  case operationType::maxMag:
154  {
155  compareScalars(mag(values), result, greaterOp<scalar>());
156  return true;
157  }
158  default:
159  {
160  return false;
161  }
162  }
163 }
164 
165 
166 template<class Type>
168 (
169  const Field<Type>& values,
170  const scalarField& weights,
171  const scalarField& V,
172  Result<Type>& result
173 ) const
174 {
175  switch (operation_)
176  {
177  case operationType::sum:
178  {
179  result.value = gSum(weights*values);
180  return true;
181  }
183  {
184  result.value = gSum(cmptMag(values));
185  return true;
186  }
188  {
189  result.value = gSum(weights*values)/max(gSum(weights), vSmall);
190  return true;
191  }
192  case operationType::volAverage:
193  {
194  result.value = gSum(weights*V*values)/max(gSum(weights*V), vSmall);
195  return true;
196  }
197  case operationType::volIntegrate:
198  {
199  result.value = gSum(weights*V*values);
200  return true;
201  }
202  case operationType::min:
203  {
204  result.value = gMin(values);
205  return true;
206  }
207  case operationType::max:
208  {
209  result.value = gMax(values);
210  return true;
211  }
212  case operationType::CoV:
213  {
214  Type meanValue = gSum(values*V)/this->V();
215 
216  const label nComp = pTraits<Type>::nComponents;
217 
218  for (direction d=0; d<nComp; ++d)
219  {
220  scalarField vals(values.component(d));
221  scalar mean = component(meanValue, d);
222  scalar& res = setComponent(result.value, d);
223 
224  res = sqrt(gSum(V*sqr(vals - mean))/this->V())/mean;
225  }
226 
227  return true;
228  }
229  case operationType::none:
230  {
231  return true;
232  }
233  default:
234  {
235  return false;
236  }
237  }
238 }
239 
240 
241 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
242 
243 template<class Type>
245 (
246  const word& fieldName,
247  const scalarField& weights,
248  const scalarField& V
249 )
250 {
251  const bool ok = validField<Type>(fieldName);
252 
253  if (ok)
254  {
255  // Get the values
256  Field<Type> values(getFieldValues<Type>(fieldName));
257 
258  // Write raw values if specified
259  if (writeFields_)
260  {
262  (
263  IOobject
264  (
265  fieldName + '_' + selectionTypeNames[selectionType()]
266  + '-' + cellSetName(),
267  time_.name(),
268  obr_,
271  ),
272  (weights*values).ref()
273  ).write();
274  }
275 
276  // Do the operation
277  if (operation_ != operationType::none)
278  {
279  // Apply scale factor
280  values *= scaleFactor_;
281 
282  bool ok = false;
283 
284  #define writeValuesFieldType(fieldType, none) \
285  ok = \
286  ok \
287  || writeValues<Type, fieldType> \
288  ( \
289  fieldName, \
290  values, \
291  weights, \
292  V \
293  );
295  #undef writeValuesFieldType
296 
297  if (!ok)
298  {
300  << "Operation " << operationTypeNames_[operation_]
301  << " not available for values of type "
303  << exit(FatalError);
304  }
305  }
306  }
307 
308  return ok;
309 }
310 
311 
312 template<class Type, class ResultType>
314 (
315  const word& fieldName,
316  const Field<Type>& values,
317  const scalarField& weights,
318  const scalarField& V
319 )
320 {
321  Result<ResultType> result({Zero, -1, -1, point::uniform(NaN)});
322 
323  if (processValues(values, weights, V, result))
324  {
325  // Add to result dictionary, over-writing any previous entry
326  resultDict_.add(fieldName, result.value, true);
327 
328  if (Pstream::master())
329  {
330  file() << tab << result.value;
331 
332  Log << " " << operationTypeNames_[operation_]
333  << "(" << cellSetName() << ") of " << fieldName
334  << " = " << result.value;
335 
336  if (result.celli != -1)
337  {
338  Log << " at location " << result.cc;
339  if (writeLocation_) file() << tab << result.cc;
340  }
341 
342  if (result.celli != -1)
343  {
344  Log << " in cell " << result.celli;
345  if (writeLocation_) file() << tab << result.celli;
346  }
347 
348  if (result.proci != -1)
349  {
350  Log << " on processor " << result.proci;
351  if (writeLocation_) file() << tab << result.proci;
352  }
353 
354  Log << endl;
355  }
356 
357  return true;
358  }
359 
360  return false;
361 }
362 
363 
364 template<class Type>
367 (
368  const Field<Type>& field
369 ) const
370 {
371  if (all())
372  {
373  return field;
374  }
375  else
376  {
377  return tmp<Field<Type>>(new Field<Type>(field, cells()));
378  }
379 }
380 
381 
382 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
Definition: Field.C:455
Generic GeometricField class.
A primitive field of type <Type> with automated input and output.
Definition: IOField.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
static Form uniform(const Cmpt &s)
Return a VectorSpace with all elements = s.
Definition: VectorSpaceI.H:168
bool writeValues(const word &fieldName, const scalarField &weights, const scalarField &V)
Templated helper function to output field values.
bool processValuesTypeType(const Field< Type > &values, const scalarField &weights, const scalarField &V, Result< Type > &result) const
Apply a Type -> Type operation to the values.
tmp< Field< Type > > getFieldValues(const word &fieldName) const
Insert field values into values list.
tmp< Field< Type > > filterField(const Field< Type > &field) const
Filter a field according to cellIds.
bool validField(const word &fieldName) const
Return true if the field name is valid.
void compareScalars(const scalarField &values, Result< scalar > &result, const Op &op) const
Apply a comparison operation to the values, returning the limiting.
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.
const fvMesh & mesh_
Reference to the fvMesh.
const objectRegistry & obr_
Reference to the region objectRegistry.
const volVectorField & C() const
Return cell centres.
bool foundObject(const word &name) const
Is the named Type in registry.
Traits class for primitives.
Definition: pTraits.H:53
virtual bool write(const bool write=true) const
Write using setting from DB.
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
const cellShapeList & cells
volScalarField & b
Definition: createFields.H:27
#define Log
Report write to Foam::Info if the local log switch is true.
tmp< VolField< Type > > average(const SurfaceField< Type > &ssf)
Area-weighted average a surfaceField creating a volField.
Definition: fvcAverage.C:46
static Type NaN()
Return a primitive with all components set to NaN.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const zero Zero
Definition: zero.H:97
Type gSum(const FieldField< Field, Type > &f)
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
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
label & setComponent(label &l, const direction)
Definition: label.H:86
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
static const char tab
Definition: Ostream.H:259
errorManip< error > abort(error &err)
Definition: errorManip.H:131
FOR_ALL_FIELD_TYPES(DefineContiguousFvWallLocationDataType)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< scalar > mag(const dimensioned< Type > &)
typename VolField< Type >::Internal VolInternalField
Definition: volFieldsFwd.H:58
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
dimensionSet cmptMag(const dimensionSet &)
Definition: dimensionSet.C:275
error FatalError
Type gMin(const FieldField< Field, Type > &f)
Type gMax(const FieldField< Field, Type > &f)
uint8_t direction
Definition: direction.H:45
dimensioned< scalar > sumMag(const DimensionedField< Type, GeoMesh > &df)
#define Op(opName, op)
Definition: ops.H:100
#define writeValuesFieldType(fieldType, none)