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  const scalar emptyVal,
83  Result<scalar>& result,
84  const Op& op
85 ) const
86 {
87  // Check if the zone on this processor is non-zero sized
88  if (values.size())
89  {
90  label i = 0;
91  forAll(values, j)
92  {
93  if (op(values[j], values[i]))
94  {
95  i = j;
96  }
97  }
98 
99  result.value = values[i];
100  result.celli = celli(i);
101  result.cc = fieldValue::mesh_.C()[result.celli];
102  }
103  else
104  {
105  result.value = emptyVal;
106  }
107 
108  result.proci = Pstream::parRun() ? Pstream::myProcNo() : -1;
109 
110  reduce
111  (
112  result,
113  [&op](const Result<scalar>& a, const Result<scalar>& b)
114  {
115  return op(a.value, b.value) ? a : b;
116  }
117  );
118 }
119 
120 
121 template<class Type, class ResultType>
123 (
124  const Field<Type>& values,
125  const scalarField& weights,
126  const scalarField& V,
127  Result<ResultType>& result
128 ) const
129 {
130  return false;
131 }
132 
133 
134 template<class Type>
136 (
137  const Field<Type>& values,
138  const scalarField& weights,
139  const scalarField& V,
140  Result<Type>& result
141 ) const
142 {
143  return processValuesTypeType(values, weights, V, result);
144 }
145 
146 
147 template<class Type>
149 (
150  const Field<Type>& values,
151  const scalarField& weights,
152  const scalarField& V,
153  Result<scalar>& result
154 ) const
155 {
156  switch (operation_)
157  {
158  case operationType::minMag:
159  {
160  compareScalars(mag(values), vGreat, result, lessOp<scalar>());
161  return true;
162  }
163  case operationType::maxMag:
164  {
165  compareScalars(mag(values), -vGreat, result, greaterOp<scalar>());
166  return true;
167  }
168  default:
169  {
170  return false;
171  }
172  }
173 }
174 
175 
176 template<class Type>
178 (
179  const Field<Type>& values,
180  const scalarField& weights,
181  const scalarField& V,
182  Result<Type>& result
183 ) const
184 {
185  switch (operation_)
186  {
187  case operationType::sum:
188  {
189  result.value = gSum(weights*values);
190  return true;
191  }
193  {
194  result.value = gSum(cmptMag(values));
195  return true;
196  }
198  {
199  result.value = gSum(weights*values)/max(gSum(weights), vSmall);
200  return true;
201  }
202  case operationType::volAverage:
203  {
204  result.value = gSum(weights*V*values)/max(gSum(weights*V), vSmall);
205  return true;
206  }
207  case operationType::volIntegrate:
208  {
209  result.value = gSum(weights*V*values);
210  return true;
211  }
212  case operationType::min:
213  {
214  result.value = gMin(values);
215  return true;
216  }
217  case operationType::max:
218  {
219  result.value = gMax(values);
220  return true;
221  }
222  case operationType::CoV:
223  {
224  Type meanValue = gSum(values*V)/this->V();
225 
226  const label nComp = pTraits<Type>::nComponents;
227 
228  for (direction d=0; d<nComp; ++d)
229  {
230  scalarField vals(values.component(d));
231  scalar mean = component(meanValue, d);
232  scalar& res = setComponent(result.value, d);
233 
234  res = sqrt(gSum(V*sqr(vals - mean))/this->V())/mean;
235  }
236 
237  return true;
238  }
239  case operationType::none:
240  {
241  return true;
242  }
243  default:
244  {
245  return false;
246  }
247  }
248 }
249 
250 
251 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
252 
253 template<class Type>
255 (
256  const word& fieldName,
257  const scalarField& weights,
258  const scalarField& V
259 )
260 {
261  const bool ok = validField<Type>(fieldName);
262 
263  if (ok)
264  {
265  // Get the values
266  Field<Type> values(getFieldValues<Type>(fieldName));
267 
268  // Write raw values if specified
269  if (writeFields_)
270  {
272  (
273  IOobject
274  (
275  fieldName + '_' + selectionTypeNames[selectionType()]
276  + '-' + cellSetName(),
277  time_.name(),
278  obr_,
281  ),
282  (weights*values).ref()
283  ).write();
284  }
285 
286  // Do the operation
287  if (operation_ != operationType::none)
288  {
289  // Apply scale factor
290  values *= scaleFactor_;
291 
292  bool ok = false;
293 
294  #define writeValuesFieldType(fieldType, none) \
295  ok = \
296  ok \
297  || writeValues<Type, fieldType> \
298  ( \
299  fieldName, \
300  values, \
301  weights, \
302  V \
303  );
305  #undef writeValuesFieldType
306 
307  if (!ok)
308  {
310  << "Operation " << operationTypeNames_[operation_]
311  << " not available for values of type "
313  << exit(FatalError);
314  }
315  }
316  }
317 
318  return ok;
319 }
320 
321 
322 template<class Type, class ResultType>
324 (
325  const word& fieldName,
326  const Field<Type>& values,
327  const scalarField& weights,
328  const scalarField& V
329 )
330 {
331  Result<ResultType> result({Zero, -1, -1, point::uniform(NaN)});
332 
333  if (processValues(values, weights, V, result))
334  {
335  // Add to result dictionary, over-writing any previous entry
336  resultDict_.add(fieldName, result.value, true);
337 
338  if (Pstream::master())
339  {
340  file() << tab << result.value;
341 
342  Log << " " << operationTypeNames_[operation_]
343  << "(" << cellSetName() << ") of " << fieldName
344  << " = " << result.value;
345 
346  if (result.celli != -1)
347  {
348  Log << " at location " << result.cc;
349  if (writeLocation_) file() << tab << result.cc;
350  }
351 
352  if (result.celli != -1)
353  {
354  Log << " in cell " << result.celli;
355  if (writeLocation_) file() << tab << result.celli;
356  }
357 
358  if (result.proci != -1)
359  {
360  Log << " on processor " << result.proci;
361  if (writeLocation_) file() << tab << result.proci;
362  }
363 
364  Log << endl;
365  }
366 
367  return true;
368  }
369 
370  return false;
371 }
372 
373 
374 template<class Type>
377 (
378  const Field<Type>& field
379 ) const
380 {
381  if (all())
382  {
383  return field;
384  }
385  else
386  {
387  return tmp<Field<Type>>(new Field<Type>(field, cells()));
388  }
389 }
390 
391 
392 // ************************************************************************* //
#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:479
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
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
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.
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.
void compareScalars(const scalarField &values, const scalar emptyVal, Result< scalar > &result, const Op &op) const
Apply a comparison operation to the values, returning the limiting.
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:334
const cellShapeList & cells
volScalarField & b
Definition: createFields.H:25
#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
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:257
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:265
errorManip< error > abort(error &err)
Definition: errorManip.H:131
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:61
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
dimensionSet cmptMag(const dimensionSet &)
Definition: dimensionSet.C:296
error FatalError
FOR_ALL_FIELD_TYPES(makeFieldSourceTypedef)
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)