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-2019 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 
42  if (obr_.foundObject<vf>(fieldName))
43  {
44  return true;
45  }
46 
47  return false;
48 }
49 
50 
51 template<class Type>
54 (
55  const word& fieldName,
56  const bool mustGet
57 ) const
58 {
60 
61  if (obr_.foundObject<vf>(fieldName))
62  {
63  return filterField(obr_.lookupObject<vf>(fieldName));
64  }
65 
66  if (mustGet)
67  {
69  << "Field " << fieldName << " not found in database"
70  << abort(FatalError);
71  }
72 
73  return tmp<Field<Type>>::New(Zero);
74 }
75 
76 
77 template<class Type>
79 (
80  const Field<Type>& values,
81  const scalarField& V,
82  const scalarField& weightField
83 ) const
84 {
85  Type result = Zero;
86  switch (operation_)
87  {
88  case opNone:
89  {
90  break;
91  }
92  case opMin:
93  {
94  result = gMin(values);
95  break;
96  }
97  case opMax:
98  {
99  result = gMax(values);
100  break;
101  }
102  case opSumMag:
103  {
104  result = gSum(cmptMag(values));
105  break;
106  }
107  case opSum:
108  case opWeightedSum:
109  {
110  if (canWeight(weightField))
111  {
112  result = gSum(weightField*values);
113  }
114  else
115  {
116  // Unweighted form
117  result = gSum(values);
118  }
119  break;
120  }
121  case opAverage:
122  case opWeightedAverage:
123  {
124  if (canWeight(weightField))
125  {
126  result =
127  gSum(weightField*values)/(gSum(weightField) + ROOTVSMALL);
128  }
129  else
130  {
131  // Unweighted form
132  const label n = returnReduce(values.size(), sumOp<label>());
133  result = gSum(values)/(scalar(n) + ROOTVSMALL);
134  }
135  break;
136  }
137  case opVolAverage:
138  case opWeightedVolAverage:
139  {
140  if (canWeight(weightField))
141  {
142  result = gSum(weightField*V*values)
143  /(gSum(weightField*V) + ROOTVSMALL);
144  }
145  else
146  {
147  // Unweighted form
148  result = gSum(V*values)/(gSum(V) + ROOTVSMALL);
149  }
150  break;
151  }
152  case opVolIntegrate:
153  case opWeightedVolIntegrate:
154  {
155  if (canWeight(weightField))
156  {
157  result = gSum(weightField*V*values);
158  }
159  else
160  {
161  // Unweighted form
162  result = gSum(V*values);
163  }
164  break;
165  }
166  case opCoV:
167  {
168  const scalar sumV = gSum(V);
169 
170  Type meanValue = gSum(V*values)/sumV;
171 
172  for (direction d=0; d < pTraits<Type>::nComponents; ++d)
173  {
174  tmp<scalarField> vals(values.component(d));
175  const scalar mean = component(meanValue, d);
176  scalar& res = setComponent(result, d);
177 
178  res = sqrt(gSum(V*sqr(vals - mean))/sumV)/(mean + ROOTVSMALL);
179  }
180 
181  break;
182  }
183  }
184 
185  return result;
186 }
187 
188 
189 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
190 
191 template<class Type>
193 (
194  const word& fieldName,
195  const scalarField& V,
196  const scalarField& weightField
197 )
198 {
199  const bool ok = validField<Type>(fieldName);
200 
201  if (ok)
202  {
203  Field<Type> values(getFieldValues<Type>(fieldName));
204 
205  if (writeFields_)
206  {
207  Field<Type> allValues(values);
208  combineFields(allValues);
209 
210  if (Pstream::master())
211  {
212  word outName = fieldName + '_' + regionTypeNames_[regionType_];
213  if (this->volRegion::regionName_ != polyMesh::defaultRegion)
214  {
215  outName = outName + '-' + this->volRegion::regionName_;
216  }
217 
219  (
220  IOobject
221  (
222  outName,
223  obr_.time().timeName(),
224  obr_,
225  IOobject::NO_READ,
226  IOobject::NO_WRITE
227  ),
228  scaleFactor_*weightField*allValues
229  ).write();
230  }
231  }
232 
233  if (operation_ != opNone)
234  {
235  // Apply scale factor
236  values *= scaleFactor_;
237 
238  Type result = processValues(values, V, weightField);
239 
240  // Write state/results information
241  const word& opName = operationTypeNames_[operation_];
242  word outName = fieldName;
243  if (this->volRegion::regionName_ != polyMesh::defaultRegion)
244  {
245  outName = this->volRegion::regionName_ + ',' + outName;
246  }
247  word resultName = opName + '(' + outName + ')';
248 
249  file()<< tab << result;
250 
251  Log << " " << opName
252  << '(' << this->volRegion::regionName_ << ") of " << fieldName
253  << " = " << result << endl;
254 
255  this->setResult(resultName, result);
256  }
257  }
258 
259  return ok;
260 }
261 
262 
263 template<class Type>
266 (
267  const Field<Type>& field
268 ) const
269 {
270  if (volRegion::vrtAll == this->volRegion::regionType())
271  {
272  return field;
273  }
274 
275  return tmp<Field<Type>>::New(field, cellIDs());
276 }
277 
278 
279 // ************************************************************************* //
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:127
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.
Definition: zero.H:128
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:79
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
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::label
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:62
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:193
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:355
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::tab
constexpr char tab
Definition: Ostream.H:371
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
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
Log
#define Log
Report write to Foam::Info if the local log switch is true.
Definition: messageStream.H:332