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-2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
34template<class Type>
36(
37 const word& fieldName
38) const
39{
41 typedef typename VolFieldType::Internal IntVolFieldType;
42
43 return
44 (
45 obr_.foundObject<VolFieldType>(fieldName)
46 || obr_.foundObject<IntVolFieldType>(fieldName)
47 );
48}
49
50
51template<class Type>
54(
55 const word& fieldName,
56 const bool mandatory
57) const
58{
60 typedef typename VolFieldType::Internal IntVolFieldType;
61
62 if (obr_.foundObject<VolFieldType>(fieldName))
63 {
64 return filterField(obr_.lookupObject<VolFieldType>(fieldName));
65 }
66 else if (obr_.foundObject<IntVolFieldType>(fieldName))
67 {
68 return filterField(obr_.lookupObject<IntVolFieldType>(fieldName));
69 }
70
71 if (mandatory)
72 {
74 << "Field " << fieldName << " not found in database" << nl
75 << abort(FatalError);
76 }
77
78 return tmp<Field<Type>>::New();
79}
80
81
82template<class Type>
84(
85 const Field<Type>& values,
86 const scalarField& V,
87 const scalarField& weightField
88) const
89{
90 Type result = Zero;
91 switch (operation_)
92 {
93 case opNone:
94 {
95 break;
96 }
97 case opMin:
98 {
99 result = gMin(values);
100 break;
101 }
102 case opMax:
103 {
104 result = gMax(values);
105 break;
106 }
107 case opSumMag:
108 {
109 result = gSum(cmptMag(values));
110 break;
111 }
112 case opSum:
113 case opWeightedSum:
114 {
115 if (is_weightedOp() && canWeight(weightField))
116 {
117 result = gSum(weightField*values);
118 }
119 else
120 {
121 // Unweighted form
122 result = gSum(values);
123 }
124 break;
125 }
126 case opAverage:
127 case opWeightedAverage:
128 {
129 if (is_weightedOp() && canWeight(weightField))
130 {
131 result =
132 gSum(weightField*values)/(gSum(weightField) + ROOTVSMALL);
133 }
134 else
135 {
136 // Unweighted form
137 const label n = returnReduce(values.size(), sumOp<label>());
138 result = gSum(values)/(scalar(n) + ROOTVSMALL);
139 }
140 break;
141 }
142 case opVolAverage:
143 case opWeightedVolAverage:
144 {
145 if (is_weightedOp() && canWeight(weightField))
146 {
147 result = gSum(weightField*V*values)
148 /(gSum(weightField*V) + ROOTVSMALL);
149 }
150 else
151 {
152 // Unweighted form
153 result = gSum(V*values)/(gSum(V) + ROOTVSMALL);
154 }
155 break;
156 }
157 case opVolIntegrate:
158 case opWeightedVolIntegrate:
159 {
160 if (is_weightedOp() && canWeight(weightField))
161 {
162 result = gSum(weightField*V*values);
163 }
164 else
165 {
166 // Unweighted form
167 result = gSum(V*values);
168 }
169 break;
170 }
171 case opCoV:
172 {
173 const scalar sumV = gSum(V);
174
175 Type meanValue = gSum(V*values)/sumV;
176
177 for (direction d=0; d < pTraits<Type>::nComponents; ++d)
178 {
179 tmp<scalarField> vals(values.component(d));
180 const scalar mean = component(meanValue, d);
181 scalar& res = setComponent(result, d);
182
183 res = sqrt(gSum(V*sqr(vals - mean))/sumV)/(mean + ROOTVSMALL);
184 }
185
186 break;
187 }
188 }
189
190 return result;
191}
192
193
194// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
195
196template<class Type>
198(
199 const word& fieldName,
200 const scalarField& V,
201 const scalarField& weightField
202)
203{
204 const bool ok = validField<Type>(fieldName);
205
206 if (ok)
207 {
208 Field<Type> values(getFieldValues<Type>(fieldName));
209
210 if (writeFields_)
211 {
212 word outName = fieldName + '_' + regionTypeNames_[regionType_];
214 {
215 outName = outName + '-' + this->volRegion::regionName_;
216 }
217
219 (
221 (
222 outName,
223 obr_.time().timeName(),
224 obr_,
227 ),
228 weightField.empty()
229 ? scaleFactor_*values
230 : scaleFactor_*weightField*values
231 ).write();
232 }
233
234 if (operation_ != opNone)
235 {
236 // Apply scale factor
237 values *= scaleFactor_;
238
239 Type result = processValues(values, V, weightField);
240
241 switch (postOperation_)
242 {
243 case postOpSqrt:
244 {
245 // sqrt: component-wise - does not change the type
246 for (direction d=0; d < pTraits<Type>::nComponents; ++d)
247 {
248 setComponent(result, d)
249 = sqrt(mag(component(result, d)));
250 }
251 break;
252 }
253 default:
254 {
255 break;
256 }
257 }
258
259 // Write state/results information
260 word prefix, suffix;
261 {
262 if (postOperation_ != postOpNone)
263 {
264 // Adjust result name to include post-operation
265 prefix += postOperationTypeNames_[postOperation_];
266 prefix += '(';
267 suffix += ')';
268 }
269
270 prefix += operationTypeNames_[operation_];
271 prefix += '(';
272 suffix += ')';
273 }
274
275 word regionPrefix;
277 {
278 regionPrefix = this->volRegion::regionName_ + ',';
279 }
280
281 word resultName = prefix + regionPrefix + fieldName + suffix;
282
283 Log << " " << prefix << this->volRegion::regionName_ << suffix
284 << " of " << fieldName << " = ";
285
286
287 // Operation or post-operation returns scalar?
288
289 scalar sresult{0};
290
291 bool alwaysScalar(operation_ & typeScalar);
292
293 if (alwaysScalar)
294 {
295 sresult = component(result, 0);
296
297 if (postOperation_ == postOpMag)
298 {
299 sresult = mag(sresult);
300 }
301 }
302 else if (postOperation_ == postOpMag)
303 {
304 sresult = mag(result);
305 alwaysScalar = true;
306 }
307
308
309 if (alwaysScalar)
310 {
311 file()<< tab << sresult;
312
313 Log << sresult << endl;
314
315 this->setResult(resultName, sresult);
316 }
317 else
318 {
319 file()<< tab << result;
320
321 Log << result << endl;
322
323 this->setResult(resultName, result);
324 }
325 }
326 }
327
328 return ok;
329}
330
331
332template<class Type>
335(
336 const Field<Type>& field
337) const
338{
339 if (this->volRegion::useAllCells())
340 {
341 return field;
342 }
343
345}
346
347
348// ************************************************************************* //
#define Log
Definition: PDRblock.C:35
label n
labelList cellIDs
Generic templated field type.
Definition: Field.H:82
Generic GeometricField class.
A primitive field of type <T> with automated input and output.
Definition: IOField.H:58
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
bool empty() const noexcept
True if the UList is empty (ie, size() is zero)
Definition: UListI.H:427
Watches for presence of the named trigger file in the case directory and signals a simulation stop (o...
Definition: abort.H:128
bool writeValues(const word &fieldName, const scalarField &V, const scalarField &weightField)
Templated helper function to output field values.
tmp< Field< Type > > getFieldValues(const word &fieldName, const bool mandatory=false) const
Insert field values into values list.
Type processValues(const Field< Type > &values, const scalarField &V, const scalarField &weightField) const
Apply the 'operation' to the values.
bool validField(const word &fieldName) const
Return true if the field name is valid.
tmp< Field< Type > > filterField(const Field< Type > &field) const
Filter a field according to cellIds.
Computes the magnitude of an input field.
Definition: mag.H:153
const objectRegistry & obr_
Reference to the region objectRegistry.
bool useAllCells() const noexcept
Use all cells, not the cellIDs.
Definition: volRegionI.H:31
wordRe regionName_
Region name (cellSet, cellZone, ...)
Definition: volRegion.H:170
bool foundObject(const word &name, const bool recursive=false) const
Is the named Type found?
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:321
virtual bool write(const bool valid=true) const
Write using setting from DB.
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
rDeltaTY field()
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Type gSum(const FieldField< Field, Type > &f)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
label & setComponent(label &val, const direction) noexcept
Non-const access to integer-type (has no components)
Definition: label.H:123
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensionedScalar sqrt(const dimensionedScalar &ds)
uint8_t direction
Definition: direction.H:56
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
error FatalError
Type gMin(const FieldField< Field, Type > &f)
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce (copy) and return value.
Type gMax(const FieldField< Field, Type > &f)
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
constexpr char tab
The tab '\t' character(0x09)
Definition: Ostream.H:52