Loading...
Searching...
No Matches
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-2023 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\*---------------------------------------------------------------------------*/
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 case typeScalar:
95 case typeWeighted:
96 case typeAbsolute:
97 {
98 break;
99 }
100 case opMin:
101 {
102 result = gMin(values);
103 break;
104 }
105 case opMax:
106 {
107 result = gMax(values);
108 break;
109 }
110 case opSumMag:
111 {
112 result = gSum(cmptMag(values));
113 break;
114 }
115 case opSum:
116 case opWeightedSum:
117 {
118 if (is_weightedOp() && canWeight(weightField))
119 {
120 result = gSum(weightField*values);
121 }
122 else
123 {
124 // Unweighted form
125 result = gSum(values);
126 }
127 break;
128 }
129 case opAverage:
130 case opWeightedAverage:
131 {
132 if (is_weightedOp() && canWeight(weightField))
133 {
134 result =
135 gSum(weightField*values)/(gSum(weightField) + ROOTVSMALL);
136 }
137 else
138 {
139 // Unweighted form
140 const label n = returnReduce(values.size(), sumOp<label>());
141 result = gSum(values)/(scalar(n) + ROOTVSMALL);
142 }
143 break;
144 }
145 case opVolAverage:
146 case opWeightedVolAverage:
147 {
148 if (is_weightedOp() && canWeight(weightField))
149 {
150 result = gSum(weightField*V*values)
151 /(gSum(weightField*V) + ROOTVSMALL);
152 }
153 else
154 {
155 // Unweighted form
156 result = gSum(V*values)/(gSum(V) + ROOTVSMALL);
157 }
158 break;
159 }
160 case opVolIntegrate:
161 case opWeightedVolIntegrate:
162 {
163 if (is_weightedOp() && canWeight(weightField))
164 {
165 result = gSum(weightField*V*values);
166 }
167 else
168 {
169 // Unweighted form
170 result = gSum(V*values);
171 }
172 break;
173 }
174 case opCoV:
175 {
176 const scalar sumV = gSum(V);
177
178 Type meanValue = gSum(V*values)/sumV;
179
180 for (direction d=0; d < pTraits<Type>::nComponents; ++d)
181 {
182 tmp<scalarField> vals(values.component(d));
183 const scalar mean = component(meanValue, d);
184 scalar& res = setComponent(result, d);
185
186 res = sqrt(gSum(V*sqr(vals - mean))/sumV)/(mean + ROOTVSMALL);
187 }
188
189 break;
190 }
191 }
192
193 return result;
194}
195
196
197// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
198
199template<class Type>
201(
202 const word& fieldName,
203 const scalarField& V,
204 const scalarField& weightField
205)
206{
207 const bool ok = validField<Type>(fieldName);
208
209 if (ok)
210 {
211 Field<Type> values(getFieldValues<Type>(fieldName));
212
213 if (writeFields_)
214 {
215 word outName = fieldName + '_' + regionTypeNames_[regionType_];
217 {
218 outName = outName + '-' + this->volRegion::regionName_;
219 }
220
222 (
224 (
225 outName,
226 obr_.time().timeName(),
227 obr_,
230 ),
231 weightField.empty()
232 ? scaleFactor_*values
233 : scaleFactor_*weightField*values
234 ).write();
235 }
236
237 if (operation_ != opNone)
238 {
239 // Apply scale factor
240 values *= scaleFactor_;
241
242 Type result = processValues(values, V, weightField);
243
244 switch (postOperation_)
245 {
246 case postOpSqrt:
247 {
248 // sqrt: component-wise - does not change the type
249 for (direction d=0; d < pTraits<Type>::nComponents; ++d)
250 {
251 setComponent(result, d)
252 = sqrt(mag(component(result, d)));
253 }
254 break;
255 }
256 default:
257 {
258 break;
259 }
260 }
261
262 // Write state/results information
263 word prefix, suffix;
264 {
265 if (postOperation_ != postOpNone)
266 {
267 // Adjust result name to include post-operation
268 prefix += postOperationTypeNames_[postOperation_];
269 prefix += '(';
270 suffix += ')';
271 }
272
273 prefix += operationTypeNames_[operation_];
274 prefix += '(';
275 suffix += ')';
276 }
277
278 word regionPrefix;
280 {
281 regionPrefix = this->volRegion::regionName_ + ',';
282 }
283
284 word resultName = prefix + regionPrefix + fieldName + suffix;
285
286 Log << " " << prefix << this->volRegion::regionName_ << suffix
287 << " of " << fieldName << " = ";
288
289
290 // Operation or post-operation returns scalar?
291
292 scalar sresult{0};
293
294 bool alwaysScalar(operation_ & typeScalar);
295
296 if (alwaysScalar)
297 {
298 sresult = component(result, 0);
299
300 if (postOperation_ == postOpMag)
301 {
302 sresult = mag(sresult);
303 }
304 }
305 else if (postOperation_ == postOpMag)
306 {
307 sresult = mag(result);
308 alwaysScalar = true;
309 }
310
311
312 if (alwaysScalar)
313 {
314 file()<< tab << sresult;
315
316 Log << sresult << endl;
317
318 this->setResult(resultName, sresult);
319 }
320 else
321 {
322 file()<< tab << result;
323
324 Log << result << endl;
325
326 this->setResult(resultName, result);
327 }
328 }
329 }
331 return ok;
332}
333
334
335template<class Type>
338(
339 const Field<Type>& field
340) const
341{
342 if (this->volRegion::useAllCells())
343 {
344 return field;
345 }
346
347 return tmp<Field<Type>>::New(field, cellIDs());
348}
349
350
351// ************************************************************************* //
#define Log
Definition PDRblock.C:28
label n
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
Generic GeometricField class.
A primitive field of type <T> with automated input and output.
Definition IOField.H:53
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
bool empty() const noexcept
True if List is empty (ie, size() is zero).
Definition UList.H:701
scalar scaleFactor_
Scaling factor.
Definition fieldValue.H:131
bool writeFields_
Flag to output field values.
Definition fieldValue.H:121
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.
postOperationType postOperation_
Optional post-evaluation operation.
Type processValues(const Field< Type > &values, const scalarField &V, const scalarField &weightField) const
Apply the 'operation' to the values.
bool canWeight(const scalarField &fld) const
True if field is non-empty on any processor.
@ postOpNone
No additional operation after calculation.
@ postOpSqrt
Component-wise sqrt after normal operation.
@ postOpMag
Component-wise mag after normal operation.
operationType operation_
Operation to apply to values.
static const Enum< operationType > operationTypeNames_
Operation type names.
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.
bool is_weightedOp() const noexcept
True if the operation variant uses a weight-field.
static const Enum< postOperationType > postOperationTypeNames_
Operation type names.
const objectRegistry & obr_
Reference to the region objectRegistry.
void setResult(const word &entryName, const Type &value)
Add result.
bool useAllCells() const noexcept
Use all cells, not the cellIDs.
Definition volRegionI.H:24
scalar V() const
Return total volume of the selected region.
Definition volRegionI.H:52
wordRe regionName_
Region name (cellSet, cellZone, ...).
Definition volRegion.H:188
static const Enum< regionTypes > regionTypeNames_
Region type names.
Definition volRegion.H:130
regionTypes regionType_
Region type.
Definition volRegion.H:183
const labelList & cellIDs() const
Return the local list of cell IDs.
Definition volRegion.C:203
virtual OFstream & file()
Return access to the file (if only 1).
Definition writeFile.C:270
static word defaultRegion
Return the default region name.
Definition polyMesh.H:406
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
rDeltaTY field()
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
List< T > values(const HashTable< T, Key, Hash > &tbl, const bool doSort=false)
List of values from HashTable, optionally sorted.
Definition HashOps.H:164
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
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:160
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
errorManip< error > abort(error &err)
Definition errorManip.H:139
uint8_t direction
Definition direction.H:49
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
Type gMin(const FieldField< Field, Type > &f)
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
Type gMax(const FieldField< Field, Type > &f)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
constexpr char tab
The tab '\t' character(0x09).
Definition Ostream.H:49