volFieldValueTemplates.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2015-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "volFieldValue.H"
30 #include "volFields.H"
31 
32 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
33 
34 template<class Type>
36 (
37  const word& fieldName
38 ) const
39 {
41  typedef typename VolFieldType::Internal IntVolFieldType;
42 
43  if
44  (
45  obr_.foundObject<VolFieldType>(fieldName)
46  || obr_.foundObject<IntVolFieldType>(fieldName)
47  )
48  {
49  return true;
50  }
51 
52  return false;
53 }
54 
55 
56 template<class Type>
59 (
60  const word& fieldName,
61  const bool mustGet
62 ) const
63 {
65  typedef typename VolFieldType::Internal IntVolFieldType;
66 
67  if (obr_.foundObject<VolFieldType>(fieldName))
68  {
69  return filterField(obr_.lookupObject<VolFieldType>(fieldName));
70  }
71  else if (obr_.foundObject<IntVolFieldType>(fieldName))
72  {
73  return filterField(obr_.lookupObject<IntVolFieldType>(fieldName));
74  }
75 
76  if (mustGet)
77  {
79  << "Field " << fieldName << " not found in database"
80  << abort(FatalError);
81  }
82 
83  return tmp<Field<Type>>::New(Zero);
84 }
85 
86 
87 template<class Type>
89 (
90  const Field<Type>& values,
91  const scalarField& V,
92  const scalarField& weightField
93 ) const
94 {
95  Type result = Zero;
96  switch (operation_)
97  {
98  case opNone:
99  {
100  break;
101  }
102  case opMin:
103  {
104  result = gMin(values);
105  break;
106  }
107  case opMax:
108  {
109  result = gMax(values);
110  break;
111  }
112  case opSumMag:
113  {
114  result = gSum(cmptMag(values));
115  break;
116  }
117  case opSum:
118  case opWeightedSum:
119  {
120  if (canWeight(weightField))
121  {
122  result = gSum(weightField*values);
123  }
124  else
125  {
126  // Unweighted form
127  result = gSum(values);
128  }
129  break;
130  }
131  case opAverage:
132  case opWeightedAverage:
133  {
134  if (canWeight(weightField))
135  {
136  result =
137  gSum(weightField*values)/(gSum(weightField) + ROOTVSMALL);
138  }
139  else
140  {
141  // Unweighted form
142  const label n = returnReduce(values.size(), sumOp<label>());
143  result = gSum(values)/(scalar(n) + ROOTVSMALL);
144  }
145  break;
146  }
147  case opVolAverage:
148  case opWeightedVolAverage:
149  {
150  if (canWeight(weightField))
151  {
152  result = gSum(weightField*V*values)
153  /(gSum(weightField*V) + ROOTVSMALL);
154  }
155  else
156  {
157  // Unweighted form
158  result = gSum(V*values)/(gSum(V) + ROOTVSMALL);
159  }
160  break;
161  }
162  case opVolIntegrate:
163  case opWeightedVolIntegrate:
164  {
165  if (canWeight(weightField))
166  {
167  result = gSum(weightField*V*values);
168  }
169  else
170  {
171  // Unweighted form
172  result = gSum(V*values);
173  }
174  break;
175  }
176  case opCoV:
177  {
178  const scalar sumV = gSum(V);
179 
180  Type meanValue = gSum(V*values)/sumV;
181 
182  for (direction d=0; d < pTraits<Type>::nComponents; ++d)
183  {
184  tmp<scalarField> vals(values.component(d));
185  const scalar mean = component(meanValue, d);
186  scalar& res = setComponent(result, d);
187 
188  res = sqrt(gSum(V*sqr(vals - mean))/sumV)/(mean + ROOTVSMALL);
189  }
190 
191  break;
192  }
193  }
194 
195  return result;
196 }
197 
198 
199 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
200 
201 template<class Type>
203 (
204  const word& fieldName,
205  const scalarField& V,
206  const scalarField& weightField
207 )
208 {
209  const bool ok = validField<Type>(fieldName);
210 
211  if (ok)
212  {
213  Field<Type> values(getFieldValues<Type>(fieldName));
214 
215  if (writeFields_)
216  {
217  Field<Type> allValues(values);
218  combineFields(allValues);
219 
220  if (Pstream::master())
221  {
222  word outName = fieldName + '_' + regionTypeNames_[regionType_];
223  if (this->volRegion::regionName_ != polyMesh::defaultRegion)
224  {
225  outName = outName + '-' + this->volRegion::regionName_;
226  }
227 
229  (
230  IOobject
231  (
232  outName,
233  obr_.time().timeName(),
234  obr_,
235  IOobject::NO_READ,
236  IOobject::NO_WRITE
237  ),
238  scaleFactor_*weightField*allValues
239  ).write();
240  }
241  }
242 
243  if (operation_ != opNone)
244  {
245  // Apply scale factor
246  values *= scaleFactor_;
247 
248  Type result = processValues(values, V, weightField);
249 
250  switch (postOperation_)
251  {
252  case postOpSqrt:
253  {
254  // sqrt: component-wise - does not change the type
255  for (direction d=0; d < pTraits<Type>::nComponents; ++d)
256  {
257  setComponent(result, d)
258  = sqrt(mag(component(result, d)));
259  }
260  break;
261  }
262  default:
263  {
264  break;
265  }
266  }
267 
268  // Write state/results information
269  word prefix, suffix;
270  {
271  if (postOperation_ != postOpNone)
272  {
273  // Adjust result name to include post-operation
274  prefix += postOperationTypeNames_[postOperation_];
275  prefix += '(';
276  suffix += ')';
277  }
278 
279  prefix += operationTypeNames_[operation_];
280  prefix += '(';
281  suffix += ')';
282  }
283 
284  word regionPrefix;
285  if (this->volRegion::regionName_ != polyMesh::defaultRegion)
286  {
287  regionPrefix = this->volRegion::regionName_ + ',';
288  }
289 
290  word resultName = prefix + regionPrefix + fieldName + suffix;
291 
292  Log << " " << prefix << this->volRegion::regionName_ << suffix
293  << " of " << fieldName << " = ";
294 
295 
296  // Operation or post-operation returns scalar?
297 
298  scalar sresult{0};
299 
300  bool alwaysScalar(operation_ & typeScalar);
301 
302  if (alwaysScalar)
303  {
304  sresult = component(result, 0);
305 
306  if (postOperation_ == postOpMag)
307  {
308  sresult = mag(sresult);
309  }
310  }
311  else if (postOperation_ == postOpMag)
312  {
313  sresult = mag(result);
314  alwaysScalar = true;
315  }
316 
317 
318  if (alwaysScalar)
319  {
320  file()<< tab << sresult;
321 
322  Log << sresult << endl;
323 
324  this->setResult(resultName, sresult);
325  }
326  else
327  {
328  file()<< tab << result;
329 
330  Log << result << endl;
331 
332  this->setResult(resultName, result);
333  }
334  }
335  }
336 
337  return ok;
338 }
339 
340 
341 template<class Type>
344 (
345  const Field<Type>& field
346 ) const
347 {
348  if (volRegion::vrtAll == this->volRegion::regionType())
349  {
350  return field;
351  }
352 
353  return tmp<Field<Type>>::New(field, cellIDs());
354 }
355 
356 
357 // ************************************************************************* //
volFields.H
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::setComponent
label & setComponent(label &l, const direction)
Definition: label.H:123
Log
#define Log
Definition: PDRblock.C:34
Foam::component
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
Definition: FieldFieldFunctions.C:44
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::IOField
A primitive field of type <T> with automated input and output.
Definition: foamVtkLagrangianWriter.H:61
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:94
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::HashTableOps::values
List< T > values(const HashTable< T, Key, Hash > &tbl, const bool doSort=false)
List of values from HashTable, optionally sorted.
Definition: HashOps.H:149
volFieldValue.H
Foam::functionObjects::fieldValues::volFieldValue::processValues
Type processValues(const Field< Type > &values, const scalarField &V, const scalarField &weightField) const
Apply the 'operation' to the values.
Definition: volFieldValueTemplates.C:89
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:594
Foam::functionObjects::fieldValues::volFieldValue::getFieldValues
tmp< Field< Type > > getFieldValues(const word &fieldName, const bool mustGet=false) const
Insert field values into values list.
Foam::sumOp
Definition: ops.H:213
Foam::cmptMag
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:400
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::functionObjects::fieldValues::volFieldValue::validField
bool validField(const word &fieldName) const
Return true if the field name is valid.
Definition: volFieldValueTemplates.C:36
Foam::Field
Generic templated field type.
Definition: Field.H:63
field
rDeltaTY field()
Foam::functionObjects::fieldValues::volFieldValue::writeValues
bool writeValues(const word &fieldName, const scalarField &V, const scalarField &weightField)
Templated helper function to output field values.
Definition: volFieldValueTemplates.C:203
Foam::functionObjects::fieldValues::volFieldValue::filterField
tmp< Field< Type > > filterField(const Field< Type > &field) const
Filter a field according to cellIds.
Foam::FatalError
error FatalError
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:137
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:372
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::tab
constexpr char tab
Definition: Ostream.H:384
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::direction
uint8_t direction
Definition: direction.H:47
Foam::gMin
Type gMin(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:593
Foam::GeometricField< Type, fvPatchField, volMesh >
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:592